24#ifndef VTK_MULTI_WRITER_HH
25#define VTK_MULTI_WRITER_HH
27#warning "This header is deprecated. Use the new vtkoutputmodule."
37#include <dune/common/fvector.hh>
38#include <dune/istl/bvector.hh>
40#include <dune/grid/io/file/vtk/vtkwriter.hh>
59template<
class Gr
idView, Dune::VTK::OutputType OutputValue = Dune::VTK::ascii >
62 enum { dim = GridView::dimension };
63 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
64 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
68 const std::string &simName =
"",
69 std::string multiFileName =
"")
71 , elementMapper_(gridView,
Dune::mcmgElementLayout())
72 , vertexMapper_(gridView,
Dune::mcmgVertexLayout())
74 simName_ = (simName.empty())?
"sim":simName;
75 multiFileName_ = multiFileName;
76 if (multiFileName_.empty()) {
77 multiFileName_ = simName_;
78 multiFileName_ +=
".pvd";
86 MPI_Comm_rank(MPI_COMM_WORLD, &commRank_);
87 MPI_Comm_size(MPI_COMM_WORLD, &commSize_);
103 {
return curWriterNum_; }
114 elementMapper_.update();
115 vertexMapper_.update();
124 if (!multiFile_.is_open()) {
125 startMultiFile_(multiFileName_);
128 curWriter_ = std::make_shared<VtkWriter>(gridView_, Dune::VTK::conforming);
130 curOutFileName_ = fileName_();
139 template <
class Scalar,
int nComp>
142 using VectorField = Dune::BlockVector<Dune::FieldVector<Scalar, nComp> >;
144 std::shared_ptr<ManagedVectorField_<VectorField> > vfs =
145 std::make_shared<ManagedVectorField_<VectorField> >(nEntities);
146 managedObjects_.push_back(vfs);
150 template <
class Scalar>
152 {
return allocateManagedBuffer<Scalar, 1>(nEntities); }
154 {
return allocateManagedBuffer<double, 1>(nEntities); }
173 template <
class DataBuffer>
176 sanitizeBuffer_(buf, nComps);
178 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
180 FunctionPtr fnPtr(
new VtkFn(name,
186 curWriter_->addVertexData(fnPtr);
204 template <
class DataBuffer>
207 sanitizeBuffer_(buf, nComps);
209 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
211 FunctionPtr fnPtr(
new VtkFn(name,
217 curWriter_->addCellData(fnPtr);
227 template <
class DataBuffer>
228 void attachDofData(DataBuffer &buf, std::string name,
bool isVertexData,
int nComps = 1)
231 attachVertexData(buf, name, nComps);
233 attachCellData(buf, name, nComps);
247 curWriter_->write(curOutFileName_.c_str(), OutputValue);
251 std::string fileName;
252 std::string suffix = fileSuffix_();
253 if (commSize_ == 1) {
254 fileName = curOutFileName_;
255 multiFile_.precision(16);
256 multiFile_ <<
" <DataSet timestep=\""
259 << fileName <<
"." << suffix <<
"\"/>\n";
261 if (commSize_ > 1 && commRank_ == 0) {
263 for (
int part=0; part < commSize_; ++part) {
264 fileName = fileName_(part);
265 multiFile_.precision(16);
266 multiFile_ <<
" <DataSet "
267 <<
" part=\"" << part <<
"\""
268 <<
" timestep=\"" << curTime_ <<
"\""
269 <<
" file=\"" << fileName <<
"." << suffix <<
"\"/>\n";
275 while (managedObjects_.begin() != managedObjects_.end()) {
276 managedObjects_.pop_front();
288 template <
class Restarter>
291 res.serializeSectionBegin(
"VTKMultiWriter");
292 res.serializeStream() << curWriterNum_ <<
"\n";
294 if (commRank_ == 0) {
297 if (multiFile_.is_open()) {
299 filePos = multiFile_.tellp();
300 multiFile_.seekp(0, std::ios::end);
301 fileLen = multiFile_.tellp();
302 multiFile_.seekp(filePos);
305 res.serializeStream() << fileLen <<
" " << filePos <<
"\n";
308 std::ifstream multiFileIn(multiFileName_.c_str());
309 char *tmp =
new char[fileLen];
310 multiFileIn.read(tmp, fileLen);
311 res.serializeStream().write(tmp, fileLen);
316 res.serializeSectionEnd();
322 template <
class Restarter>
325 res.deserializeSectionBegin(
"VTKMultiWriter");
326 res.deserializeStream() >> curWriterNum_;
328 if (commRank_ == 0) {
330 std::getline(res.deserializeStream(), dummy);
333 size_t filePos, fileLen;
334 res.deserializeStream() >> fileLen >> filePos;
335 std::getline(res.deserializeStream(), dummy);
336 if (multiFile_.is_open())
340 multiFile_.open(multiFileName_.c_str());
342 char *tmp =
new char[fileLen];
343 res.deserializeStream().read(tmp, fileLen);
344 multiFile_.write(tmp, fileLen);
348 multiFile_.seekp(filePos);
352 std::getline(res.deserializeStream(), tmp);
354 res.deserializeSectionEnd();
358 std::string fileName_()
360 std::ostringstream oss;
361 oss << simName_ <<
"-" << std::setw(5) << std::setfill(
'0') << curWriterNum_;
365 std::string fileName_(
int rank)
368 std::ostringstream oss;
369 oss <<
"s" << std::setw(4) << std::setfill(
'0') << commSize_
370 <<
"-p" << std::setw(4) << std::setfill(
'0') << rank
371 <<
"-" << simName_ <<
"-"
372 << std::setw(5) << curWriterNum_;
378 std::ostringstream oss;
379 oss << simName_ <<
"-" << std::setw(5) << std::setfill(
'0') << curWriterNum_;
384 std::string fileSuffix_()
386 return (GridView::dimension == 1)?
"vtp":
"vtu";
390 void startMultiFile_(
const std::string &multiFileName)
393 if (commRank_ == 0) {
395 multiFile_.open(multiFileName.c_str());
397 "<?xml version=\"1.0\"?>\n"
398 "<VTKFile type=\"Collection\"\n"
400 " byte_order=\"LittleEndian\"\n"
401 " compressor=\"vtkZLibDataCompressor\">\n"
406 void finishMultiFile_()
409 if (commRank_ == 0) {
411 std::ofstream::pos_type pos = multiFile_.tellp();
415 multiFile_.seekp(pos);
422 template <
class DataBuffer>
423 void sanitizeBuffer_(DataBuffer &b,
int nComps)
425 for (
unsigned int i = 0; i < b.size(); ++i) {
426 for (
int j = 0; j < nComps; ++j) {
431 if (std::abs(b[i][j])
432 < std::numeric_limits<float>::min())
459 virtual ~ManagedObject_()
464 class ManagedVectorField_ :
public ManagedObject_
467 ManagedVectorField_(
int size)
475 const GridView gridView_;
476 ElementMapper elementMapper_;
477 VertexMapper vertexMapper_;
479 std::string simName_;
480 std::ofstream multiFile_;
481 std::string multiFileName_;
486 std::shared_ptr<VtkWriter> curWriter_;
488 std::string curOutFileName_;
491 std::list<std::shared_ptr<ManagedObject_> > managedObjects_;
Some templates to wrap the valgrind macros.
Provides a vector-valued function using Dune::FieldVectors as elements.
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
Definition: common/pdesolver.hh:35
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:61
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:289
int curWriterNum() const
Returns the number of the current VTK file.
Definition: vtkmultiwriter.hh:102
Dune::BlockVector< Dune::FieldVector< Scalar, nComp > > * allocateManagedBuffer(int nEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:140
void beginWrite(double t)
Called when ever a new time step or a new grid must be written.
Definition: vtkmultiwriter.hh:122
Dune::BlockVector< Dune::FieldVector< double, 1 > > * allocateManagedBuffer(int nEntities)
Definition: vtkmultiwriter.hh:153
void attachVertexData(DataBuffer &buf, std::string name, int nComps=1)
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:174
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition: vtkmultiwriter.hh:112
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:323
Dune::VTKWriter< GridView > VtkWriter
Definition: vtkmultiwriter.hh:66
void attachDofData(DataBuffer &buf, std::string name, bool isVertexData, int nComps=1)
Add data associated with degrees of freedom to the output.
Definition: vtkmultiwriter.hh:228
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > * allocateManagedBuffer(int nEntities)
Definition: vtkmultiwriter.hh:151
void attachCellData(DataBuffer &buf, std::string name, int nComps=1)
Add a cell centered quantity to the output.
Definition: vtkmultiwriter.hh:205
~VtkMultiWriter()
Definition: vtkmultiwriter.hh:91
VtkMultiWriter(const GridView &gridView, const std::string &simName="", std::string multiFileName="")
Definition: vtkmultiwriter.hh:67
void endWrite(bool onlyDiscard=false)
Finalizes the current writer.
Definition: vtkmultiwriter.hh:243
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtknestedfunction.hh:43
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:313