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>
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};
92 multiWriter_.multiFile_.precision(16);
93 multiWriter_.multiFile_ <<
" <DataSet timestep=\"" << multiWriter_.curTime_ <<
"\" file=\""
101 enum { dim = GridView::dimension };
103 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
104 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
107 using Scalar = BaseOutputWriter::Scalar;
108 using Vector = BaseOutputWriter::Vector;
109 using Tensor = BaseOutputWriter::Tensor;
110 using ScalarBuffer = BaseOutputWriter::ScalarBuffer;
111 using VectorBuffer = BaseOutputWriter::VectorBuffer;
112 using TensorBuffer = BaseOutputWriter::TensorBuffer;
114 using VtkWriter = Dune::VTKWriter<GridView>;
115 using FunctionPtr = std::shared_ptr< Dune::VTKFunction< GridView > >;
118 const GridView& gridView,
119 const std::string& outputDir,
120 const std::string&
simName =
"",
122 : gridView_(gridView)
123 , elementMapper_(gridView, Dune::mcmgElementLayout())
124 , vertexMapper_(gridView, Dune::mcmgVertexLayout())
125 , curWriter_(
nullptr)
129 outputDir_ = outputDir;
135 if (multiFileName_.empty())
136 multiFileName_ = outputDir_+
"/"+simName_+
".pvd";
138 commRank_ = gridView.comm().rank();
139 commSize_ = gridView.comm().size();
156 {
return curWriterNum_; }
167 elementMapper_.update(gridView_);
168 vertexMapper_.update(gridView_);
176 if (!multiFile_.is_open()) {
177 startMultiFile_(multiFileName_);
186 curOutFileName_ = fileName_();
188 curWriter_ =
new VtkWriter(gridView_, Dune::VTK::conforming);
201 managedScalarBuffers_.push_back(
buf);
214 for (
size_t i = 0; i <
numOuter; ++ i)
217 managedVectorBuffers_.push_back(
buf);
239 sanitizeScalarBuffer_(
buf);
247 curWriter_->addVertexData(
fnPtr);
267 sanitizeScalarBuffer_(
buf);
275 curWriter_->addCellData(
fnPtr);
296 sanitizeVectorBuffer_(
buf);
304 curWriter_->addVertexData(
fnPtr);
315 std::ostringstream
oss;
324 curWriter_->addVertexData(
fnPtr);
345 sanitizeVectorBuffer_(
buf);
353 curWriter_->addCellData(
fnPtr);
364 std::ostringstream
oss;
373 curWriter_->addCellData(
fnPtr);
387 auto tasklet = std::make_shared<WriteDataTasklet>(*
this);
402 template <
class Restarter>
405 res.serializeSectionBegin(
"VTKMultiWriter");
406 res.serializeStream() << curWriterNum_ <<
"\n";
408 if (commRank_ == 0) {
411 if (multiFile_.is_open()) {
414 multiFile_.seekp(0, std::ios::end);
430 res.serializeSectionEnd();
436 template <
class Restarter>
439 res.deserializeSectionBegin(
"VTKMultiWriter");
440 res.deserializeStream() >> curWriterNum_;
442 if (commRank_ == 0) {
444 std::getline(
res.deserializeStream(),
dummy);
450 std::getline(
res.deserializeStream(),
dummy);
451 if (multiFile_.is_open())
455 multiFile_.open(multiFileName_.c_str());
459 multiFile_.write(tmp,
fileLen);
467 std::getline(
res.deserializeStream(), tmp);
469 res.deserializeSectionEnd();
473 std::string fileName_()
476 std::ostringstream
oss;
477 oss << simName_ <<
"-"
478 << std::setw(5) << std::setfill(
'0') << curWriterNum_;
482 std::string fileSuffix_()
483 {
return (GridView::dimension == 1) ?
"vtp" :
"vtu"; }
488 if (commRank_ == 0) {
491 multiFile_ <<
"<?xml version=\"1.0\"?>\n"
492 "<VTKFile type=\"Collection\"\n"
494 " byte_order=\"LittleEndian\"\n"
495 " compressor=\"vtkZLibDataCompressor\">\n"
500 void finishMultiFile_()
503 if (commRank_ == 0) {
505 std::ofstream::pos_type pos = multiFile_.tellp();
506 multiFile_ <<
" </Collection>\n"
508 multiFile_.seekp(pos);
515 void sanitizeScalarBuffer_(ScalarBuffer&)
520 void sanitizeVectorBuffer_(VectorBuffer&)
526 void releaseBuffers_()
530 curWriter_ =
nullptr;
531 while (managedScalarBuffers_.begin() != managedScalarBuffers_.end()) {
532 delete managedScalarBuffers_.front();
533 managedScalarBuffers_.pop_front();
535 while (managedVectorBuffers_.begin() != managedVectorBuffers_.end()) {
536 delete managedVectorBuffers_.front();
537 managedVectorBuffers_.pop_front();
541 const GridView gridView_;
542 ElementMapper elementMapper_;
543 VertexMapper vertexMapper_;
545 std::string outputDir_;
546 std::string simName_;
547 std::ofstream multiFile_;
548 std::string multiFileName_;
553 VtkWriter *curWriter_;
555 std::string curOutFileName_;
558 std::list<ScalarBuffer *> managedScalarBuffers_;
559 std::list<VectorBuffer *> managedVectorBuffers_;
561 TaskletRunner taskletRunner_;
The base class for all output writers.
The base class for all output writers.
Definition baseoutputwriter.hh:44
The base class for tasklets.
Definition tasklets.hpp:44
void barrier()
Make sure that all tasklets have been completed after this method has been called.
Definition tasklets.cpp:134
void dispatch(std::shared_ptr< TaskletInterface > tasklet)
Add a new tasklet.
Definition tasklets.cpp:101
Simplifies writing multi-file VTK datasets.
Definition vtkmultiwriter.hh:66
void attachVectorElementData(VectorBuffer &buf, std::string name) override
Add a element centered quantity to the output.
Definition vtkmultiwriter.hh:343
void endWrite(bool onlyDiscard=false) override
Finalizes the current writer.
Definition vtkmultiwriter.hh:384
void attachTensorVertexData(TensorBuffer &buf, std::string name) override
Add a finished vertex-centered tensor field to the output.
Definition vtkmultiwriter.hh:310
void attachTensorElementData(TensorBuffer &buf, std::string name) override
Add a finished element-centered tensor field to the output.
Definition vtkmultiwriter.hh:359
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition vtkmultiwriter.hh:437
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:211
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition vtkmultiwriter.hh:165
void beginWrite(double t) override
Called whenever a new time step must be written.
Definition vtkmultiwriter.hh:174
ScalarBuffer * allocateManagedScalarBuffer(size_t numEntities)
Allocate a managed buffer for a scalar field.
Definition vtkmultiwriter.hh:198
void attachScalarElementData(ScalarBuffer &buf, std::string name) override
Add a element centered quantity to the output.
Definition vtkmultiwriter.hh:265
void attachScalarVertexData(ScalarBuffer &buf, std::string name) override
Add a finished vertex centered vector field to the output.
Definition vtkmultiwriter.hh:237
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition vtkmultiwriter.hh:403
void attachVectorVertexData(VectorBuffer &buf, std::string name) override
Add a finished vertex centered vector field to the output.
Definition vtkmultiwriter.hh:294
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
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
Provides a mechanism to dispatch work to separate threads.
Provides a vector-valued function using Dune::FieldVectors as elements.
Provides a tensor-valued function using Dune::FieldMatrix objects as elements.
Provides a vector-valued function using Dune::FieldVectors as elements.