28#ifndef EWOMS_VTK_MULTI_WRITER_HH
29#define EWOMS_VTK_MULTI_WRITER_HH
38#include <opm/material/common/Valgrind.hpp>
40#include <dune/common/fvector.hh>
41#include <dune/common/version.hh>
42#include <dune/istl/bvector.hh>
43#include <dune/grid/io/file/vtk/vtkwriter.hh>
64template <
class Gr
idView,
int vtkFormat>
71 : multiWriter_(multiWriter)
78 if (multiWriter_.commSize_ > 1)
79 fileName = multiWriter_.curWriter_->pwrite(multiWriter_.curOutFileName_,
80 multiWriter_.outputDir_,
82 static_cast<Dune::VTK::OutputType
>(vtkFormat));
84 fileName = multiWriter_.curWriter_->write(multiWriter_.outputDir_ +
"/" + multiWriter_.curOutFileName_,
85 static_cast<Dune::VTK::OutputType
>(vtkFormat));
90 const std::filesystem::path fullPath{fileName};
91 const std::string localFileName = fullPath.filename();
92 multiWriter_.multiFile_.precision(16);
93 multiWriter_.multiFile_ <<
" <DataSet timestep=\"" << multiWriter_.curTime_ <<
"\" file=\""
94 << localFileName <<
"\"/>\n";
101 enum { dim = GridView::dimension };
103 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
104 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
115 using FunctionPtr = std::shared_ptr< Dune::VTKFunction< GridView > >;
118 const GridView& gridView,
119 const std::string& outputDir,
120 const std::string& simName =
"",
121 std::string multiFileName =
"")
122 : gridView_(gridView)
123 , elementMapper_(gridView,
Dune::mcmgElementLayout())
124 , vertexMapper_(gridView,
Dune::mcmgVertexLayout())
125 , curWriter_(nullptr)
127 , taskletRunner_(asyncWriting?1:0)
129 outputDir_ = outputDir;
133 simName_ = (simName.empty()) ?
"sim" : simName;
134 multiFileName_ = multiFileName;
135 if (multiFileName_.empty())
136 multiFileName_ = outputDir_+
"/"+simName_+
".pvd";
138 commRank_ = gridView.comm().rank();
139 commSize_ = gridView.comm().size();
156 {
return curWriterNum_; }
167#if DUNE_VERSION_NEWER(DUNE_GRID, 2, 8)
168 elementMapper_.update(gridView_);
169 vertexMapper_.update(gridView_);
171 elementMapper_.update();
172 vertexMapper_.update();
181 if (!multiFile_.is_open()) {
182 startMultiFile_(multiFileName_);
191 curOutFileName_ = fileName_();
193 curWriter_ =
new VtkWriter(gridView_, Dune::VTK::conforming);
206 managedScalarBuffers_.push_back(buf);
219 for (
size_t i = 0; i < numOuter; ++ i)
220 (*buf)[i].resize(numInner);
222 managedVectorBuffers_.push_back(buf);
244 sanitizeScalarBuffer_(buf);
252 curWriter_->addVertexData(fnPtr);
272 sanitizeScalarBuffer_(buf);
280 curWriter_->addCellData(fnPtr);
301 sanitizeVectorBuffer_(buf);
309 curWriter_->addVertexData(fnPtr);
319 for (
unsigned colIdx = 0; colIdx < buf[0].N(); ++colIdx) {
320 std::ostringstream oss;
321 oss << name <<
"[" << colIdx <<
"]";
329 curWriter_->addVertexData(fnPtr);
350 sanitizeVectorBuffer_(buf);
358 curWriter_->addCellData(fnPtr);
368 for (
unsigned colIdx = 0; colIdx < buf[0].N(); ++colIdx) {
369 std::ostringstream oss;
370 oss << name <<
"[" << colIdx <<
"]";
378 curWriter_->addCellData(fnPtr);
392 auto tasklet = std::make_shared<WriteDataTasklet>(*
this);
407 template <
class Restarter>
410 res.serializeSectionBegin(
"VTKMultiWriter");
411 res.serializeStream() << curWriterNum_ <<
"\n";
413 if (commRank_ == 0) {
414 std::streamsize fileLen = 0;
415 std::streamoff filePos = 0;
416 if (multiFile_.is_open()) {
418 filePos = multiFile_.tellp();
419 multiFile_.seekp(0, std::ios::end);
420 fileLen = multiFile_.tellp();
421 multiFile_.seekp(filePos);
424 res.serializeStream() << fileLen <<
" " << filePos <<
"\n";
427 std::ifstream multiFileIn(multiFileName_.c_str());
428 char *tmp =
new char[fileLen];
429 multiFileIn.read(tmp,
static_cast<long>(fileLen));
430 res.serializeStream().write(tmp, fileLen);
435 res.serializeSectionEnd();
441 template <
class Restarter>
444 res.deserializeSectionBegin(
"VTKMultiWriter");
445 res.deserializeStream() >> curWriterNum_;
447 if (commRank_ == 0) {
449 std::getline(res.deserializeStream(), dummy);
452 std::streamoff filePos;
453 std::streamsize fileLen;
454 res.deserializeStream() >> fileLen >> filePos;
455 std::getline(res.deserializeStream(), dummy);
456 if (multiFile_.is_open())
460 multiFile_.open(multiFileName_.c_str());
462 char *tmp =
new char[fileLen];
463 res.deserializeStream().read(tmp, fileLen);
464 multiFile_.write(tmp, fileLen);
468 multiFile_.seekp(filePos);
472 std::getline(res.deserializeStream(), tmp);
474 res.deserializeSectionEnd();
478 std::string fileName_()
481 std::ostringstream oss;
482 oss << simName_ <<
"-"
483 << std::setw(5) << std::setfill(
'0') << curWriterNum_;
487 std::string fileSuffix_()
488 {
return (GridView::dimension == 1) ?
"vtp" :
"vtu"; }
490 void startMultiFile_(
const std::string& multiFileName)
493 if (commRank_ == 0) {
495 multiFile_.open(multiFileName.c_str());
496 multiFile_ <<
"<?xml version=\"1.0\"?>\n"
497 "<VTKFile type=\"Collection\"\n"
499 " byte_order=\"LittleEndian\"\n"
500 " compressor=\"vtkZLibDataCompressor\">\n"
505 void finishMultiFile_()
508 if (commRank_ == 0) {
510 std::ofstream::pos_type pos = multiFile_.tellp();
511 multiFile_ <<
" </Collection>\n"
513 multiFile_.seekp(pos);
531 void releaseBuffers_()
535 curWriter_ =
nullptr;
536 while (managedScalarBuffers_.begin() != managedScalarBuffers_.end()) {
537 delete managedScalarBuffers_.front();
538 managedScalarBuffers_.pop_front();
540 while (managedVectorBuffers_.begin() != managedVectorBuffers_.end()) {
541 delete managedVectorBuffers_.front();
542 managedVectorBuffers_.pop_front();
546 const GridView gridView_;
547 ElementMapper elementMapper_;
548 VertexMapper vertexMapper_;
550 std::string outputDir_;
551 std::string simName_;
552 std::ofstream multiFile_;
553 std::string multiFileName_;
560 std::string curOutFileName_;
563 std::list<ScalarBuffer *> managedScalarBuffers_;
564 std::list<VectorBuffer *> managedVectorBuffers_;
566 TaskletRunner taskletRunner_;
The base class for all output writers.
Definition: baseoutputwriter.hh:44
Dune::DynamicVector< double > Vector
Definition: baseoutputwriter.hh:47
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:51
Dune::DynamicMatrix< double > Tensor
Definition: baseoutputwriter.hh:48
std::vector< Vector > VectorBuffer
Definition: baseoutputwriter.hh:50
std::vector< Scalar > ScalarBuffer
Definition: baseoutputwriter.hh:49
double Scalar
Definition: baseoutputwriter.hh:46
The base class for tasklets.
Definition: tasklets.hh:47
void barrier()
Make sure that all tasklets have been completed after this method has been called.
Definition: tasklets.hh:278
void dispatch(std::shared_ptr< TaskletInterface > tasklet)
Add a new tasklet.
Definition: tasklets.hh:230
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
void endWrite(bool onlyDiscard=false)
Finalizes the current writer.
Definition: vtkmultiwriter.hh:389
void attachTensorVertexData(TensorBuffer &buf, std::string name)
Add a finished vertex-centered tensor field to the output.
Definition: vtkmultiwriter.hh:315
VtkMultiWriter(bool asyncWriting, const GridView &gridView, const std::string &outputDir, const std::string &simName="", std::string multiFileName="")
Definition: vtkmultiwriter.hh:117
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:442
void attachTensorElementData(TensorBuffer &buf, std::string name)
Add a finished element-centered tensor field to the output.
Definition: vtkmultiwriter.hh:364
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: vtkmultiwriter.hh:110
int curWriterNum() const
Returns the number of the current VTK file.
Definition: vtkmultiwriter.hh:155
VectorBuffer * allocateManagedVectorBuffer(size_t numOuter, size_t numInner)
Allocate a managed buffer for a vector field.
Definition: vtkmultiwriter.hh:216
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition: vtkmultiwriter.hh:165
void beginWrite(double t)
Called whenever a new time step must be written.
Definition: vtkmultiwriter.hh:179
~VtkMultiWriter()
Definition: vtkmultiwriter.hh:142
ScalarBuffer * allocateManagedScalarBuffer(size_t numEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:203
void attachScalarVertexData(ScalarBuffer &buf, std::string name)
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:242
void attachVectorElementData(VectorBuffer &buf, std::string name)
Add a element centered quantity to the output.
Definition: vtkmultiwriter.hh:348
BaseOutputWriter::VectorBuffer VectorBuffer
Definition: vtkmultiwriter.hh:111
void attachScalarElementData(ScalarBuffer &buf, std::string name)
Add a element centered quantity to the output.
Definition: vtkmultiwriter.hh:270
Dune::VTKWriter< GridView > VtkWriter
Definition: vtkmultiwriter.hh:114
std::shared_ptr< Dune::VTKFunction< GridView > > FunctionPtr
Definition: vtkmultiwriter.hh:115
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:408
void attachVectorVertexData(VectorBuffer &buf, std::string name)
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:299
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:47
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtkvectorfunction.hh:50
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:37
Provides a mechanism to dispatch work to separate threads.