24#ifndef VTK_MULTI_WRITER_HH
25#define VTK_MULTI_WRITER_HH
27#warning "This header is deprecated. Use the new vtkoutputmodule."
35#include <dune/common/fvector.hh>
36#include <dune/istl/bvector.hh>
38#include <dune/grid/io/file/vtk/vtkwriter.hh>
39#include <dune/grid/io/file/vtk/function.hh>
55template <
class Gr
idView,
class Mapper,
class Buffer>
58 enum { dim = GridView::dimension };
59 using ctype =
typename GridView::ctype;
60 using Element =
typename GridView::template Codim<0>::Entity;
64 const GridView &gridView,
76 assert(std::size_t(buf_.size()) == std::size_t(mapper_.size()));
79 virtual std::string
name ()
const
86 const Element &element,
87 const Dune::FieldVector< ctype, dim > &xi)
const
92 idx = mapper_.index(element);
94 else if (codim_ == dim) {
99 int n = element.subEntities(dim);
101 for (
int i=0; i < n; ++i)
103 Dune::FieldVector<ctype,dim> local = referenceElement(element).position(i,dim);
105 if (local.infinity_norm()<min)
107 min = local.infinity_norm();
113 idx = mapper_.subIndex(element, imin, codim_);
116 DUNE_THROW(Dune::InvalidStateException,
117 "Only element and vertex based vector "
118 " fields are supported so far.");
120 double val = buf_[idx][mycomp];
122 if (abs(val) < std::numeric_limits<float>::min())
129 const std::string name_;
130 const GridView gridView_;
131 const Mapper &mapper_;
146template<
class Gr
idView, Dune::VTK::OutputType OutputValue = Dune::VTK::ascii >
149 enum { dim = GridView::dimension };
150 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
151 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
155 const std::string &simName =
"",
156 std::string multiFileName =
"")
157 : gridView_(gridView)
158 , elementMapper_(gridView,
Dune::mcmgElementLayout())
159 , vertexMapper_(gridView,
Dune::mcmgVertexLayout())
161 simName_ = (simName.empty())?
"sim":simName;
162 multiFileName_ = multiFileName;
163 if (multiFileName_.empty()) {
164 multiFileName_ = simName_;
165 multiFileName_ +=
".pvd";
173 MPI_Comm_rank(MPI_COMM_WORLD, &commRank_);
174 MPI_Comm_size(MPI_COMM_WORLD, &commSize_);
190 {
return curWriterNum_; }
201 if constexpr (Deprecated::hasUpdateGridView<ElementMapper, GridView>())
202 elementMapper_.update(gridView_);
204 Deprecated::update(elementMapper_);
206 if constexpr (Deprecated::hasUpdateGridView<VertexMapper, GridView>())
207 vertexMapper_.update(gridView_);
209 Deprecated::update(vertexMapper_);
218 if (!multiFile_.is_open()) {
219 startMultiFile_(multiFileName_);
222 curWriter_ = std::make_shared<VtkWriter>(gridView_, Dune::VTK::conforming);
224 curOutFileName_ = fileName_();
233 template <
class Scalar,
int nComp>
236 using VectorField = Dune::BlockVector<Dune::FieldVector<Scalar, nComp> >;
238 std::shared_ptr<ManagedVectorField_<VectorField> > vfs =
239 std::make_shared<ManagedVectorField_<VectorField> >(nEntities);
240 managedObjects_.push_back(vfs);
244 template <
class Scalar>
246 {
return allocateManagedBuffer<Scalar, 1>(nEntities); }
248 {
return allocateManagedBuffer<double, 1>(nEntities); }
267 template <
class DataBuffer>
270 sanitizeBuffer_(buf, nComps);
272 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
274 FunctionPtr fnPtr(
new VtkFn(name,
280 curWriter_->addVertexData(fnPtr);
298 template <
class DataBuffer>
301 sanitizeBuffer_(buf, nComps);
303 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
305 FunctionPtr fnPtr(
new VtkFn(name,
311 curWriter_->addCellData(fnPtr);
321 template <
class DataBuffer>
322 void attachDofData(DataBuffer &buf, std::string name,
bool isVertexData,
int nComps = 1)
325 attachVertexData(buf, name, nComps);
327 attachCellData(buf, name, nComps);
341 curWriter_->write(curOutFileName_.c_str(), OutputValue);
345 std::string fileName;
346 std::string suffix = fileSuffix_();
347 if (commSize_ == 1) {
348 fileName = curOutFileName_;
349 multiFile_.precision(16);
350 multiFile_ <<
" <DataSet timestep=\""
353 << fileName <<
"." << suffix <<
"\"/>\n";
355 if (commSize_ > 1 && commRank_ == 0) {
357 for (
int part=0; part < commSize_; ++part) {
358 fileName = fileName_(part);
359 multiFile_.precision(16);
360 multiFile_ <<
" <DataSet "
361 <<
" part=\"" << part <<
"\""
362 <<
" timestep=\"" << curTime_ <<
"\""
363 <<
" file=\"" << fileName <<
"." << suffix <<
"\"/>\n";
369 while (managedObjects_.begin() != managedObjects_.end()) {
370 managedObjects_.pop_front();
382 template <
class Restarter>
385 res.serializeSectionBegin(
"VTKMultiWriter");
386 res.serializeStream() << curWriterNum_ <<
"\n";
388 if (commRank_ == 0) {
391 if (multiFile_.is_open()) {
393 filePos = multiFile_.tellp();
394 multiFile_.seekp(0, std::ios::end);
395 fileLen = multiFile_.tellp();
396 multiFile_.seekp(filePos);
399 res.serializeStream() << fileLen <<
" " << filePos <<
"\n";
402 std::ifstream multiFileIn(multiFileName_.c_str());
403 char *tmp =
new char[fileLen];
404 multiFileIn.read(tmp, fileLen);
405 res.serializeStream().write(tmp, fileLen);
410 res.serializeSectionEnd();
416 template <
class Restarter>
419 res.deserializeSectionBegin(
"VTKMultiWriter");
420 res.deserializeStream() >> curWriterNum_;
422 if (commRank_ == 0) {
424 std::getline(res.deserializeStream(), dummy);
427 size_t filePos, fileLen;
428 res.deserializeStream() >> fileLen >> filePos;
429 std::getline(res.deserializeStream(), dummy);
430 if (multiFile_.is_open())
434 multiFile_.open(multiFileName_.c_str());
436 char *tmp =
new char[fileLen];
437 res.deserializeStream().read(tmp, fileLen);
438 multiFile_.write(tmp, fileLen);
442 multiFile_.seekp(filePos);
446 std::getline(res.deserializeStream(), tmp);
448 res.deserializeSectionEnd();
452 std::string fileName_()
454 std::ostringstream oss;
455 oss << simName_ <<
"-" << std::setw(5) << std::setfill(
'0') << curWriterNum_;
459 std::string fileName_(
int rank)
462 std::ostringstream oss;
463 oss <<
"s" << std::setw(4) << std::setfill(
'0') << commSize_
464 <<
"-p" << std::setw(4) << std::setfill(
'0') << rank
465 <<
"-" << simName_ <<
"-"
466 << std::setw(5) << curWriterNum_;
472 std::ostringstream oss;
473 oss << simName_ <<
"-" << std::setw(5) << std::setfill(
'0') << curWriterNum_;
478 std::string fileSuffix_()
480 return (GridView::dimension == 1)?
"vtp":
"vtu";
484 void startMultiFile_(
const std::string &multiFileName)
487 if (commRank_ == 0) {
489 multiFile_.open(multiFileName.c_str());
491 "<?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();
509 multiFile_.seekp(pos);
515 template <
class DataBuffer>
516 void sanitizeBuffer_(DataBuffer &b,
int nComps)
518 for (
unsigned int i = 0; i < b.size(); ++i) {
519 for (
int j = 0; j < nComps; ++j) {
522 if (std::abs(b[i][j])
523 < std::numeric_limits<float>::min())
550 virtual ~ManagedObject_()
555 class ManagedVectorField_ :
public ManagedObject_
558 ManagedVectorField_(
int size)
566 const GridView gridView_;
567 ElementMapper elementMapper_;
568 VertexMapper vertexMapper_;
570 std::string simName_;
571 std::ofstream multiFile_;
572 std::string multiFileName_;
577 std::shared_ptr<VtkWriter> curWriter_;
579 std::string curOutFileName_;
582 std::list<std::shared_ptr<ManagedObject_> > managedObjects_;
Definition: common/pdesolver.hh:36
Provides a vector-valued function using Dune::FieldVectors as elements. DEPRECATED will be removed on...
Definition: vtkmultiwriter.hh:57
virtual int ncomps() const
Definition: vtkmultiwriter.hh:82
virtual double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtkmultiwriter.hh:85
virtual std::string name() const
Definition: vtkmultiwriter.hh:79
VtkNestedFunction(std::string name, const GridView &gridView, const Mapper &mapper, const Buffer &buf, int codim, int numComp)
Definition: vtkmultiwriter.hh:63
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:148
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:383
int curWriterNum() const
Returns the number of the current VTK file.
Definition: vtkmultiwriter.hh:189
Dune::BlockVector< Dune::FieldVector< Scalar, nComp > > * allocateManagedBuffer(int nEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:234
void beginWrite(double t)
Called when ever a new time step or a new grid must be written.
Definition: vtkmultiwriter.hh:216
Dune::BlockVector< Dune::FieldVector< double, 1 > > * allocateManagedBuffer(int nEntities)
Definition: vtkmultiwriter.hh:247
void attachVertexData(DataBuffer &buf, std::string name, int nComps=1)
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:268
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition: vtkmultiwriter.hh:199
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:417
Dune::VTKWriter< GridView > VtkWriter
Definition: vtkmultiwriter.hh:153
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:322
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > * allocateManagedBuffer(int nEntities)
Definition: vtkmultiwriter.hh:245
void attachCellData(DataBuffer &buf, std::string name, int nComps=1)
Add a cell centered quantity to the output.
Definition: vtkmultiwriter.hh:299
~VtkMultiWriter()
Definition: vtkmultiwriter.hh:178
VtkMultiWriter(const GridView &gridView, const std::string &simName="", std::string multiFileName="")
Definition: vtkmultiwriter.hh:154
void endWrite(bool onlyDiscard=false)
Finalizes the current writer.
Definition: vtkmultiwriter.hh:337
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:307