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>
53template <
class Gr
idView,
class Mapper,
class Buffer>
56 enum { dim = GridView::dimension };
57 using ctype =
typename GridView::ctype;
58 using Element =
typename GridView::template Codim<0>::Entity;
62 const GridView &gridView,
74 assert(std::size_t(buf_.size()) == std::size_t(mapper_.size()));
77 virtual std::string
name ()
const
84 const Element &element,
85 const Dune::FieldVector< ctype, dim > &xi)
const
90 idx = mapper_.index(element);
92 else if (codim_ == dim) {
97 int n = element.subEntities(dim);
99 for (
int i=0; i < n; ++i)
101 Dune::FieldVector<ctype,dim> local = referenceElement(element).position(i,dim);
103 if (local.infinity_norm()<min)
105 min = local.infinity_norm();
111 idx = mapper_.subIndex(element, imin, codim_);
114 DUNE_THROW(Dune::InvalidStateException,
115 "Only element and vertex based vector "
116 " fields are supported so far.");
118 double val = buf_[idx][mycomp];
120 if (abs(val) < std::numeric_limits<float>::min())
127 const std::string name_;
128 const GridView gridView_;
129 const Mapper &mapper_;
144template<
class Gr
idView, Dune::VTK::OutputType OutputValue = Dune::VTK::ascii >
147 enum { dim = GridView::dimension };
148 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
149 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
153 const std::string &simName =
"",
154 std::string multiFileName =
"")
155 : gridView_(gridView)
156 , elementMapper_(gridView,
Dune::mcmgElementLayout())
157 , vertexMapper_(gridView,
Dune::mcmgVertexLayout())
159 simName_ = (simName.empty())?
"sim":simName;
160 multiFileName_ = multiFileName;
161 if (multiFileName_.empty()) {
162 multiFileName_ = simName_;
163 multiFileName_ +=
".pvd";
171 MPI_Comm_rank(MPI_COMM_WORLD, &commRank_);
172 MPI_Comm_size(MPI_COMM_WORLD, &commSize_);
188 {
return curWriterNum_; }
199 elementMapper_.update();
200 vertexMapper_.update();
209 if (!multiFile_.is_open()) {
210 startMultiFile_(multiFileName_);
213 curWriter_ = std::make_shared<VtkWriter>(gridView_, Dune::VTK::conforming);
215 curOutFileName_ = fileName_();
224 template <
class Scalar,
int nComp>
227 using VectorField = Dune::BlockVector<Dune::FieldVector<Scalar, nComp> >;
229 std::shared_ptr<ManagedVectorField_<VectorField> > vfs =
230 std::make_shared<ManagedVectorField_<VectorField> >(nEntities);
231 managedObjects_.push_back(vfs);
235 template <
class Scalar>
237 {
return allocateManagedBuffer<Scalar, 1>(nEntities); }
239 {
return allocateManagedBuffer<double, 1>(nEntities); }
258 template <
class DataBuffer>
261 sanitizeBuffer_(buf, nComps);
263 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
265 FunctionPtr fnPtr(
new VtkFn(name,
271 curWriter_->addVertexData(fnPtr);
289 template <
class DataBuffer>
292 sanitizeBuffer_(buf, nComps);
294 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
296 FunctionPtr fnPtr(
new VtkFn(name,
302 curWriter_->addCellData(fnPtr);
312 template <
class DataBuffer>
313 void attachDofData(DataBuffer &buf, std::string name,
bool isVertexData,
int nComps = 1)
316 attachVertexData(buf, name, nComps);
318 attachCellData(buf, name, nComps);
332 curWriter_->write(curOutFileName_.c_str(), OutputValue);
336 std::string fileName;
337 std::string suffix = fileSuffix_();
338 if (commSize_ == 1) {
339 fileName = curOutFileName_;
340 multiFile_.precision(16);
341 multiFile_ <<
" <DataSet timestep=\""
344 << fileName <<
"." << suffix <<
"\"/>\n";
346 if (commSize_ > 1 && commRank_ == 0) {
348 for (
int part=0; part < commSize_; ++part) {
349 fileName = fileName_(part);
350 multiFile_.precision(16);
351 multiFile_ <<
" <DataSet "
352 <<
" part=\"" << part <<
"\""
353 <<
" timestep=\"" << curTime_ <<
"\""
354 <<
" file=\"" << fileName <<
"." << suffix <<
"\"/>\n";
360 while (managedObjects_.begin() != managedObjects_.end()) {
361 managedObjects_.pop_front();
373 template <
class Restarter>
376 res.serializeSectionBegin(
"VTKMultiWriter");
377 res.serializeStream() << curWriterNum_ <<
"\n";
379 if (commRank_ == 0) {
382 if (multiFile_.is_open()) {
384 filePos = multiFile_.tellp();
385 multiFile_.seekp(0, std::ios::end);
386 fileLen = multiFile_.tellp();
387 multiFile_.seekp(filePos);
390 res.serializeStream() << fileLen <<
" " << filePos <<
"\n";
393 std::ifstream multiFileIn(multiFileName_.c_str());
394 char *tmp =
new char[fileLen];
395 multiFileIn.read(tmp, fileLen);
396 res.serializeStream().write(tmp, fileLen);
401 res.serializeSectionEnd();
407 template <
class Restarter>
410 res.deserializeSectionBegin(
"VTKMultiWriter");
411 res.deserializeStream() >> curWriterNum_;
413 if (commRank_ == 0) {
415 std::getline(res.deserializeStream(), dummy);
418 size_t filePos, fileLen;
419 res.deserializeStream() >> fileLen >> filePos;
420 std::getline(res.deserializeStream(), dummy);
421 if (multiFile_.is_open())
425 multiFile_.open(multiFileName_.c_str());
427 char *tmp =
new char[fileLen];
428 res.deserializeStream().read(tmp, fileLen);
429 multiFile_.write(tmp, fileLen);
433 multiFile_.seekp(filePos);
437 std::getline(res.deserializeStream(), tmp);
439 res.deserializeSectionEnd();
443 std::string fileName_()
445 std::ostringstream oss;
446 oss << simName_ <<
"-" << std::setw(5) << std::setfill(
'0') << curWriterNum_;
450 std::string fileName_(
int rank)
453 std::ostringstream oss;
454 oss <<
"s" << std::setw(4) << std::setfill(
'0') << commSize_
455 <<
"-p" << std::setw(4) << std::setfill(
'0') << rank
456 <<
"-" << simName_ <<
"-"
457 << std::setw(5) << curWriterNum_;
463 std::ostringstream oss;
464 oss << simName_ <<
"-" << std::setw(5) << std::setfill(
'0') << curWriterNum_;
469 std::string fileSuffix_()
471 return (GridView::dimension == 1)?
"vtp":
"vtu";
475 void startMultiFile_(
const std::string &multiFileName)
478 if (commRank_ == 0) {
480 multiFile_.open(multiFileName.c_str());
482 "<?xml version=\"1.0\"?>\n"
483 "<VTKFile type=\"Collection\"\n"
485 " byte_order=\"LittleEndian\"\n"
486 " compressor=\"vtkZLibDataCompressor\">\n"
491 void finishMultiFile_()
494 if (commRank_ == 0) {
496 std::ofstream::pos_type pos = multiFile_.tellp();
500 multiFile_.seekp(pos);
506 template <
class DataBuffer>
507 void sanitizeBuffer_(DataBuffer &b,
int nComps)
509 for (
unsigned int i = 0; i < b.size(); ++i) {
510 for (
int j = 0; j < nComps; ++j) {
513 if (std::abs(b[i][j])
514 < std::numeric_limits<float>::min())
541 virtual ~ManagedObject_()
546 class ManagedVectorField_ :
public ManagedObject_
549 ManagedVectorField_(
int size)
557 const GridView gridView_;
558 ElementMapper elementMapper_;
559 VertexMapper vertexMapper_;
561 std::string simName_;
562 std::ofstream multiFile_;
563 std::string multiFileName_;
568 std::shared_ptr<VtkWriter> curWriter_;
570 std::string curOutFileName_;
573 std::list<std::shared_ptr<ManagedObject_> > managedObjects_;
Definition: common/pdesolver.hh:35
Provides a vector-valued function using Dune::FieldVectors as elements. DEPRECATED will be removed on...
Definition: vtkmultiwriter.hh:55
virtual int ncomps() const
Definition: vtkmultiwriter.hh:80
virtual double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtkmultiwriter.hh:83
virtual std::string name() const
Definition: vtkmultiwriter.hh:77
VtkNestedFunction(std::string name, const GridView &gridView, const Mapper &mapper, const Buffer &buf, int codim, int numComp)
Definition: vtkmultiwriter.hh:61
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:146
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:374
int curWriterNum() const
Returns the number of the current VTK file.
Definition: vtkmultiwriter.hh:187
Dune::BlockVector< Dune::FieldVector< Scalar, nComp > > * allocateManagedBuffer(int nEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:225
void beginWrite(double t)
Called when ever a new time step or a new grid must be written.
Definition: vtkmultiwriter.hh:207
Dune::BlockVector< Dune::FieldVector< double, 1 > > * allocateManagedBuffer(int nEntities)
Definition: vtkmultiwriter.hh:238
void attachVertexData(DataBuffer &buf, std::string name, int nComps=1)
Add a finished vertex centered vector field to the output.
Definition: vtkmultiwriter.hh:259
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition: vtkmultiwriter.hh:197
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:408
Dune::VTKWriter< GridView > VtkWriter
Definition: vtkmultiwriter.hh:151
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:313
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > * allocateManagedBuffer(int nEntities)
Definition: vtkmultiwriter.hh:236
void attachCellData(DataBuffer &buf, std::string name, int nComps=1)
Add a cell centered quantity to the output.
Definition: vtkmultiwriter.hh:290
~VtkMultiWriter()
Definition: vtkmultiwriter.hh:176
VtkMultiWriter(const GridView &gridView, const std::string &simName="", std::string multiFileName="")
Definition: vtkmultiwriter.hh:152
void endWrite(bool onlyDiscard=false)
Finalizes the current writer.
Definition: vtkmultiwriter.hh:328
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:305