vtkmultiwriter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_VTK_MULTI_WRITER_HH
29#define EWOMS_VTK_MULTI_WRITER_HH
30
31#include <dune/common/fvector.hh>
32#include <dune/common/version.hh>
33#include <dune/grid/io/file/vtk/vtkwriter.hh>
34#include <dune/istl/bvector.hh>
35
36#include <opm/material/common/Valgrind.hpp>
37
42
44
45#include <cstddef>
46#include <filesystem>
47#include <fstream>
48#include <memory>
49#include <sstream>
50#include <string>
51#include <string_view>
52#include <vector>
53
54namespace Opm {
55
63template <class GridView, Dune::VTK::OutputType vtkFormat>
65{
66 class WriteDataTasklet : public TaskletInterface
67 {
68 public:
69 explicit WriteDataTasklet(VtkMultiWriter& multiWriter)
70 : multiWriter_(multiWriter)
71 {}
72
73 void run() override final
74 {
75 std::string fileName;
76 // write the actual data as vtu or vtp (plus the pieces file in the parallel case)
77 if (multiWriter_.commSize_ > 1) {
78 fileName = multiWriter_.curWriter_->pwrite(/*name=*/multiWriter_.curOutFileName_,
79 /*path=*/multiWriter_.outputDir_,
80 /*extendPath=*/"",
81 vtkFormat);
82 }
83 else {
84 fileName = multiWriter_.curWriter_->write(/*name=*/multiWriter_.outputDir_ + "/" +
85 multiWriter_.curOutFileName_,
86 vtkFormat);
87 }
88
89 // determine name to write into the multi-file for the
90 // current time step
91 // The file names in the pvd file are relative, the path should therefore be stripped.
92 const std::filesystem::path fullPath{fileName};
93 const std::string localFileName = fullPath.filename();
94 multiWriter_.multiFile_.precision(16);
95 multiWriter_.multiFile_ << " <DataSet timestep=\"" << multiWriter_.curTime_ << "\" file=\""
96 << localFileName << "\"/>\n";
97 }
98
99 private:
100 VtkMultiWriter& multiWriter_;
101 };
102
103 enum { dim = GridView::dimension };
104
105 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
106 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
107
108public:
115
116 using VtkWriter = Dune::VTKWriter<GridView>;
117
118 VtkMultiWriter(bool asyncWriting,
119 const GridView& gridView,
120 const std::string& outputDir,
121 const std::string& simName = "",
122 const std::string& multiFileName = "")
123 : gridView_(gridView)
124 , elementMapper_(gridView, Dune::mcmgElementLayout())
125 , vertexMapper_(gridView, Dune::mcmgVertexLayout())
126 , outputDir_(outputDir.empty() ? "." : outputDir)
127 , simName_(simName.empty() ? "sim" : simName)
128 , multiFileName_(multiFileName.empty() ? outputDir_ + "/" + simName_ + ".pvd" : multiFileName)
129 , commSize_(gridView.comm().size())
130 , commRank_(gridView.comm().rank())
131 , curWriterNum_(0)
132 , taskletRunner_(/*numThreads=*/asyncWriting ? 1 : 0)
133 {}
134
136 {
137 taskletRunner_.barrier();
138 releaseBuffers_();
139 finishMultiFile_();
140
141 if (commRank_ == 0) {
142 multiFile_.close();
143 }
144 }
145
149 int curWriterNum() const
150 { return curWriterNum_; }
151
160 {
161 elementMapper_.update(gridView_);
162 vertexMapper_.update(gridView_);
163 }
164
168 void beginWrite(double t) override
169 {
170 if (!multiFile_.is_open()) {
171 startMultiFile_(multiFileName_);
172 }
173
174 // make sure that all previous output has been written and no other thread
175 // accesses the memory used as the target for the extracted quantities
176 taskletRunner_.barrier();
177 releaseBuffers_();
178
179 curTime_ = t;
180 curOutFileName_ = fileName_();
181
182 curWriter_ = std::make_unique<VtkWriter>(gridView_, Dune::VTK::conforming);
183 ++curWriterNum_;
184 }
185
193 {
194 managedScalarBuffers_.emplace_back(std::make_unique<ScalarBuffer>(numEntities));
195 return managedScalarBuffers_.back().get();
196 }
197
204 VectorBuffer* allocateManagedVectorBuffer(std::size_t numOuter, std::size_t numInner)
205 {
206 managedVectorBuffers_.emplace_back(std::make_unique<VectorBuffer>(numOuter));
207 for (auto& buffer : *managedVectorBuffers_.back()) {
208 buffer.resize(numInner);
209 }
210
211 return managedVectorBuffers_.back().get();
212 }
213
230 void attachScalarVertexData(ScalarBuffer& buf, std::string_view name) override
231 {
232 sanitizeScalarBuffer_(buf);
233
235 curWriter_->addVertexData(std::make_shared<VtkFn>(name,
236 gridView_,
237 vertexMapper_,
238 buf,
239 /*codim=*/dim));
240 }
241
257 void attachScalarElementData(ScalarBuffer& buf, std::string_view name) override
258 {
259 sanitizeScalarBuffer_(buf);
260
262 curWriter_->addCellData(std::make_shared<VtkFn>(name,
263 gridView_,
264 elementMapper_,
265 buf,
266 /*codim=*/0));
267 }
268
285 void attachVectorVertexData(VectorBuffer& buf, std::string_view name) override
286 {
287 sanitizeVectorBuffer_(buf);
288
290 curWriter_->addVertexData(std::make_shared<VtkFn>(name,
291 gridView_,
292 vertexMapper_,
293 buf,
294 /*codim=*/dim));
295 }
296
300 void attachTensorVertexData(TensorBuffer& buf, std::string_view name) override
301 {
303
304 for (unsigned colIdx = 0; colIdx < buf[0].N(); ++colIdx) {
305 std::ostringstream oss;
306 oss << name << "[" << colIdx << "]";
307 curWriter_->addVertexData(std::make_shared<VtkFn>(oss.str(),
308 gridView_,
309 vertexMapper_,
310 buf,
311 /*codim=*/dim,
312 colIdx));
313 }
314 }
315
331 void attachVectorElementData(VectorBuffer& buf, std::string_view name) override
332 {
333 sanitizeVectorBuffer_(buf);
334
336 curWriter_->addCellData(std::make_shared<VtkFn>(name,
337 gridView_,
338 elementMapper_,
339 buf,
340 /*codim=*/0));
341 }
342
346 void attachTensorElementData(TensorBuffer& buf, std::string_view name) override
347 {
349
350 for (unsigned colIdx = 0; colIdx < buf[0].N(); ++colIdx) {
351 std::ostringstream oss;
352 oss << name << "[" << colIdx << "]";
353 curWriter_->addCellData(std::make_shared<VtkFn>(oss.str(),
354 gridView_,
355 elementMapper_,
356 buf,
357 /*codim=*/0,
358 colIdx));
359 }
360 }
361
369 void endWrite(bool onlyDiscard) override
370 {
371 if (!onlyDiscard) {
372 auto tasklet = std::make_shared<WriteDataTasklet>(*this);
373 taskletRunner_.dispatch(tasklet);
374 }
375 else {
376 --curWriterNum_;
377 }
378
379 // temporarily write the closing XML mumbo-jumbo to the mashup
380 // file so that the data set can be loaded even if the
381 // simulation is aborted (or not yet finished)
382 finishMultiFile_();
383 }
384
388 template <class Restarter>
389 void serialize(Restarter& res)
390 {
391 res.serializeSectionBegin("VTKMultiWriter");
392 res.serializeStream() << curWriterNum_ << "\n";
393
394 if (commRank_ == 0) {
395 std::streamsize fileLen = 0;
396 std::streamoff filePos = 0;
397 if (multiFile_.is_open()) {
398 // write the meta file into the restart file
399 filePos = multiFile_.tellp();
400 multiFile_.seekp(0, std::ios::end);
401 fileLen = multiFile_.tellp();
402 multiFile_.seekp(filePos);
403 }
404
405 res.serializeStream() << fileLen << " " << filePos << "\n";
406
407 if (fileLen > 0) {
408 std::ifstream multiFileIn(multiFileName_.c_str());
409 std::vector<char> tmp(fileLen);
410 multiFileIn.read(tmp.data(), static_cast<long>(fileLen));
411 res.serializeStream().write(tmp.data(), fileLen);
412 }
413 }
414
415 res.serializeSectionEnd();
416 }
417
421 template <class Restarter>
422 void deserialize(Restarter& res)
423 {
424 res.deserializeSectionBegin("VTKMultiWriter");
425 res.deserializeStream() >> curWriterNum_;
426
427 if (commRank_ == 0) {
428 std::string dummy;
429 std::getline(res.deserializeStream(), dummy);
430
431 // recreate the meta file from the restart file
432 std::streamoff filePos;
433 std::streamsize fileLen;
434 res.deserializeStream() >> fileLen >> filePos;
435 std::getline(res.deserializeStream(), dummy);
436 if (multiFile_.is_open()) {
437 multiFile_.close();
438 }
439
440 if (fileLen > 0) {
441 multiFile_.open(multiFileName_);
442
443 std::vector<char> tmp(fileLen);
444 res.deserializeStream().read(tmp.data(), fileLen);
445 multiFile_.write(tmp.data(), fileLen);
446 }
447
448 multiFile_.seekp(filePos);
449 }
450 else {
451 std::string tmp;
452 std::getline(res.deserializeStream(), tmp);
453 }
454 res.deserializeSectionEnd();
455 }
456
457private:
458 std::string fileName_() const
459 {
460 // use a new file name for each time step
461 std::ostringstream oss;
462 oss << simName_ << "-"
463 << std::setw(5) << std::setfill('0') << curWriterNum_;
464 return oss.str();
465 }
466
467 static std::string fileSuffix_()
468 { return (GridView::dimension == 1) ? "vtp" : "vtu"; }
469
470 void startMultiFile_(const std::string& multiFileName)
471 {
472 // only the first process writes to the multi-file
473 if (commRank_ == 0) {
474 // generate one meta vtk-file holding the individual time steps
475 multiFile_.open(multiFileName);
476 multiFile_ << "<?xml version=\"1.0\"?>\n"
477 "<VTKFile type=\"Collection\"\n"
478 " version=\"0.1\"\n"
479 " byte_order=\"LittleEndian\"\n"
480 " compressor=\"vtkZLibDataCompressor\">\n"
481 " <Collection>\n";
482 }
483 }
484
485 void finishMultiFile_()
486 {
487 // only the first process writes to the multi-file
488 if (commRank_ == 0) {
489 // make sure that we always have a working meta file
490 std::ofstream::pos_type pos = multiFile_.tellp();
491 multiFile_ << " </Collection>\n"
492 "</VTKFile>\n";
493 multiFile_.seekp(pos);
494 multiFile_.flush();
495 }
496 }
497
498 // make sure the field is well defined if running under valgrind
499 // and make sure that all values can be displayed by paraview
500 void sanitizeScalarBuffer_(ScalarBuffer&)
501 {
502 // nothing to do: this is done by VtkScalarFunction
503 }
504
505 void sanitizeVectorBuffer_(VectorBuffer&)
506 {
507 // nothing to do: this is done by VtkVectorFunction
508 }
509
510 // release the memory occupied by all buffer objects managed by the multi-writer
511 void releaseBuffers_()
512 {
513 // discard managed objects and the current VTK writer
514 curWriter_.reset();
515 managedScalarBuffers_.clear();
516 managedVectorBuffers_.clear();
517 }
518
519 const GridView gridView_;
520 ElementMapper elementMapper_;
521 VertexMapper vertexMapper_;
522
523 const std::string outputDir_;
524 const std::string simName_;
525 std::ofstream multiFile_;
526 const std::string multiFileName_;
527
528 const int commSize_; // number of processes in the communicator
529 const int commRank_; // rank of the current process in the communicator
530
531 std::unique_ptr<VtkWriter> curWriter_;
532 double curTime_{};
533 std::string curOutFileName_;
534 int curWriterNum_;
535
536 std::vector<std::unique_ptr<ScalarBuffer>> managedScalarBuffers_;
537 std::vector<std::unique_ptr<VectorBuffer>> managedVectorBuffers_;
538
539 TaskletRunner taskletRunner_;
540};
541
542} // namespace Opm
543
544#endif
The base class for all output writers.
Definition: baseoutputwriter.hh:46
Dune::DynamicVector< double > Vector
Definition: baseoutputwriter.hh:50
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:52
Dune::DynamicMatrix< double > Tensor
Definition: baseoutputwriter.hh:49
std::vector< Vector > VectorBuffer
Definition: baseoutputwriter.hh:53
std::vector< Scalar > ScalarBuffer
Definition: baseoutputwriter.hh:51
double Scalar
Definition: baseoutputwriter.hh:48
The base class for tasklets.
Definition: tasklets.hpp:45
void barrier()
Make sure that all tasklets have been completed after this method has been called.
void dispatch(std::shared_ptr< TaskletInterface > tasklet)
Add a new tasklet.
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:65
VectorBuffer * allocateManagedVectorBuffer(std::size_t numOuter, std::size_t numInner)
Allocate a managed buffer for a vector field.
Definition: vtkmultiwriter.hh:204
void attachTensorElementData(TensorBuffer &buf, std::string_view name) override
Add a finished element-centered tensor field to the output.
Definition: vtkmultiwriter.hh:346
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:422
void endWrite(bool onlyDiscard) override
Finalizes the current writer.
Definition: vtkmultiwriter.hh:369
void attachScalarVertexData(ScalarBuffer &buf, std::string_view name) override
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:230
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: vtkmultiwriter.hh:112
int curWriterNum() const
Returns the number of the current VTK file.
Definition: vtkmultiwriter.hh:149
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition: vtkmultiwriter.hh:159
void beginWrite(double t) override
Called whenever a new time step must be written.
Definition: vtkmultiwriter.hh:168
void attachVectorVertexData(VectorBuffer &buf, std::string_view name) override
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:285
ScalarBuffer * allocateManagedScalarBuffer(std::size_t numEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:192
void attachVectorElementData(VectorBuffer &buf, std::string_view name) override
Add a element centered quantity to the output.
Definition: vtkmultiwriter.hh:331
BaseOutputWriter::VectorBuffer VectorBuffer
Definition: vtkmultiwriter.hh:114
Dune::VTKWriter< GridView > VtkWriter
Definition: vtkmultiwriter.hh:116
void attachScalarElementData(ScalarBuffer &buf, std::string_view name) override
Add a element centered quantity to the output.
Definition: vtkmultiwriter.hh:257
VtkMultiWriter(bool asyncWriting, const GridView &gridView, const std::string &outputDir, const std::string &simName="", const std::string &multiFileName="")
Definition: vtkmultiwriter.hh:118
void attachTensorVertexData(TensorBuffer &buf, std::string_view name) override
Add a finished vertex-centered tensor field to the output.
Definition: vtkmultiwriter.hh:300
~VtkMultiWriter() override
Definition: vtkmultiwriter.hh:135
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:389
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtkscalarfunction.hh:50
Provides a tensor-valued function using Dune::FieldMatrix objects as elements.
Definition: vtktensorfunction.hh:49
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtkvectorfunction.hh:48
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:39
Provides a mechanism to dispatch work to separate threads.