28#ifndef EWOMS_VTK_MULTI_WRITER_HH
29#define EWOMS_VTK_MULTI_WRITER_HH
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>
36#include <opm/material/common/Valgrind.hpp>
63template <
class Gr
idView, Dune::VTK::OutputType vtkFormat>
70 : multiWriter_(multiWriter)
73 void run()
override final
77 if (multiWriter_.commSize_ > 1) {
78 fileName = multiWriter_.curWriter_->pwrite(multiWriter_.curOutFileName_,
79 multiWriter_.outputDir_,
84 fileName = multiWriter_.curWriter_->write(multiWriter_.outputDir_ +
"/" +
85 multiWriter_.curOutFileName_,
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";
103 enum { dim = GridView::dimension };
105 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
106 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
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())
132 , taskletRunner_(asyncWriting ? 1 : 0)
141 if (commRank_ == 0) {
150 {
return curWriterNum_; }
161 elementMapper_.update(gridView_);
162 vertexMapper_.update(gridView_);
170 if (!multiFile_.is_open()) {
171 startMultiFile_(multiFileName_);
180 curOutFileName_ = fileName_();
182 curWriter_ = std::make_unique<VtkWriter>(gridView_, Dune::VTK::conforming);
194 managedScalarBuffers_.emplace_back(std::make_unique<ScalarBuffer>(numEntities));
195 return managedScalarBuffers_.back().get();
206 managedVectorBuffers_.emplace_back(std::make_unique<VectorBuffer>(numOuter));
207 for (
auto& buffer : *managedVectorBuffers_.back()) {
208 buffer.resize(numInner);
211 return managedVectorBuffers_.back().get();
232 sanitizeScalarBuffer_(buf);
235 curWriter_->addVertexData(std::make_shared<VtkFn>(name,
259 sanitizeScalarBuffer_(buf);
262 curWriter_->addCellData(std::make_shared<VtkFn>(name,
287 sanitizeVectorBuffer_(buf);
290 curWriter_->addVertexData(std::make_shared<VtkFn>(name,
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(),
333 sanitizeVectorBuffer_(buf);
336 curWriter_->addCellData(std::make_shared<VtkFn>(name,
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(),
372 auto tasklet = std::make_shared<WriteDataTasklet>(*
this);
388 template <
class Restarter>
391 res.serializeSectionBegin(
"VTKMultiWriter");
392 res.serializeStream() << curWriterNum_ <<
"\n";
394 if (commRank_ == 0) {
395 std::streamsize fileLen = 0;
396 std::streamoff filePos = 0;
397 if (multiFile_.is_open()) {
399 filePos = multiFile_.tellp();
400 multiFile_.seekp(0, std::ios::end);
401 fileLen = multiFile_.tellp();
402 multiFile_.seekp(filePos);
405 res.serializeStream() << fileLen <<
" " << filePos <<
"\n";
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);
415 res.serializeSectionEnd();
421 template <
class Restarter>
424 res.deserializeSectionBegin(
"VTKMultiWriter");
425 res.deserializeStream() >> curWriterNum_;
427 if (commRank_ == 0) {
429 std::getline(res.deserializeStream(), dummy);
432 std::streamoff filePos;
433 std::streamsize fileLen;
434 res.deserializeStream() >> fileLen >> filePos;
435 std::getline(res.deserializeStream(), dummy);
436 if (multiFile_.is_open()) {
441 multiFile_.open(multiFileName_);
443 std::vector<char> tmp(fileLen);
444 res.deserializeStream().read(tmp.data(), fileLen);
445 multiFile_.write(tmp.data(), fileLen);
448 multiFile_.seekp(filePos);
452 std::getline(res.deserializeStream(), tmp);
454 res.deserializeSectionEnd();
458 std::string fileName_()
const
461 std::ostringstream oss;
462 oss << simName_ <<
"-"
463 << std::setw(5) << std::setfill(
'0') << curWriterNum_;
467 static std::string fileSuffix_()
468 {
return (GridView::dimension == 1) ?
"vtp" :
"vtu"; }
470 void startMultiFile_(
const std::string& multiFileName)
473 if (commRank_ == 0) {
475 multiFile_.open(multiFileName);
476 multiFile_ <<
"<?xml version=\"1.0\"?>\n"
477 "<VTKFile type=\"Collection\"\n"
479 " byte_order=\"LittleEndian\"\n"
480 " compressor=\"vtkZLibDataCompressor\">\n"
485 void finishMultiFile_()
488 if (commRank_ == 0) {
490 std::ofstream::pos_type pos = multiFile_.tellp();
491 multiFile_ <<
" </Collection>\n"
493 multiFile_.seekp(pos);
511 void releaseBuffers_()
515 managedScalarBuffers_.clear();
516 managedVectorBuffers_.clear();
519 const GridView gridView_;
520 ElementMapper elementMapper_;
521 VertexMapper vertexMapper_;
523 const std::string outputDir_;
524 const std::string simName_;
525 std::ofstream multiFile_;
526 const std::string multiFileName_;
531 std::unique_ptr<VtkWriter> curWriter_;
533 std::string curOutFileName_;
536 std::vector<std::unique_ptr<ScalarBuffer>> managedScalarBuffers_;
537 std::vector<std::unique_ptr<VectorBuffer>> managedVectorBuffers_;
539 TaskletRunner taskletRunner_;
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.