3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
vtkmultiwriter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef VTK_MULTI_WRITER_HH
25#define VTK_MULTI_WRITER_HH
26
27#warning "This header is deprecated. Use the new vtkoutputmodule."
28
29#include <iostream>
30#include <limits>
31#include <list>
32#include <memory>
33#include <string>
34
35#include <dune/common/fvector.hh>
36#include <dune/istl/bvector.hh>
37
38#include <dune/grid/io/file/vtk/vtkwriter.hh>
39#include <dune/grid/io/file/vtk/function.hh>
40
42
43#if HAVE_MPI
44#include <mpi.h>
45#endif
46
47namespace Dumux {
48
55template <class GridView, class Mapper, class Buffer>
56class VtkNestedFunction : public Dune::VTKFunction<GridView>
57{
58 enum { dim = GridView::dimension };
59 using ctype = typename GridView::ctype;
60 using Element = typename GridView::template Codim<0>::Entity;
61
62public:
64 const GridView &gridView,
65 const Mapper &mapper,
66 const Buffer &buf,
67 int codim,
68 int numComp)
69 : name_(name)
70 , gridView_(gridView)
71 , mapper_(mapper)
72 , buf_(buf)
73 , codim_(codim)
74 , numComp_(numComp)
75 {
76 assert(std::size_t(buf_.size()) == std::size_t(mapper_.size()));
77 }
78
79 virtual std::string name () const
80 { return name_; }
81
82 virtual int ncomps() const
83 { return numComp_; }
84
85 virtual double evaluate(int mycomp,
86 const Element &element,
87 const Dune::FieldVector< ctype, dim > &xi) const
88 {
89 int idx;
90 if (codim_ == 0) {
91 // cells. map element to the index
92 idx = mapper_.index(element);
93 }
94 else if (codim_ == dim) {
95 // find vertex which is closest to xi in local
96 // coordinates. This code is based on Dune::P1VTKFunction
97 double min=1e100;
98 int imin=-1;
99 int n = element.subEntities(dim);
100
101 for (int i=0; i < n; ++i)
102 {
103 Dune::FieldVector<ctype,dim> local = referenceElement(element).position(i,dim);
104 local -= xi;
105 if (local.infinity_norm()<min)
106 {
107 min = local.infinity_norm();
108 imin = i;
109 }
110 }
111
112 // map vertex to an index
113 idx = mapper_.subIndex(element, imin, codim_);
114 }
115 else
116 DUNE_THROW(Dune::InvalidStateException,
117 "Only element and vertex based vector "
118 " fields are supported so far.");
119
120 double val = buf_[idx][mycomp];
121 using std::abs;
122 if (abs(val) < std::numeric_limits<float>::min())
123 val = 0;
124
125 return val;
126 }
127
128private:
129 const std::string name_;
130 const GridView gridView_;
131 const Mapper &mapper_;
132 const Buffer &buf_;
133 int codim_;
134 int numComp_;
135};
136
146template<class GridView, Dune::VTK::OutputType OutputValue = Dune::VTK::ascii >
147class [[deprecated("Use VtkOutputModule instead!")]] VtkMultiWriter
148{
149 enum { dim = GridView::dimension };
150 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
151 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
152public:
153 using VtkWriter = Dune::VTKWriter<GridView>;
154 VtkMultiWriter(const GridView &gridView,
155 const std::string &simName = "",
156 std::string multiFileName = "")
157 : gridView_(gridView)
158 , elementMapper_(gridView, Dune::mcmgElementLayout())
159 , vertexMapper_(gridView, Dune::mcmgVertexLayout())
160 {
161 simName_ = (simName.empty())?"sim":simName;
162 multiFileName_ = multiFileName;
163 if (multiFileName_.empty()) {
164 multiFileName_ = simName_;
165 multiFileName_ += ".pvd";
166 }
167
168 curWriterNum_ = 0;
169
170 commRank_ = 0;
171 commSize_ = 1;
172#if HAVE_MPI
173 MPI_Comm_rank(MPI_COMM_WORLD, &commRank_);
174 MPI_Comm_size(MPI_COMM_WORLD, &commSize_);
175#endif
176 }
177
179 {
180 finishMultiFile_();
181
182 if (commRank_ == 0)
183 multiFile_.close();
184 }
185
189 int curWriterNum() const
190 { return curWriterNum_; }
191
200 {
201 if constexpr (Deprecated::hasUpdateGridView<ElementMapper, GridView>())
202 elementMapper_.update(gridView_);
203 else
204 Deprecated::update(elementMapper_);
205
206 if constexpr (Deprecated::hasUpdateGridView<VertexMapper, GridView>())
207 vertexMapper_.update(gridView_);
208 else
209 Deprecated::update(vertexMapper_);
210 }
211
216 void beginWrite(double t)
217 {
218 if (!multiFile_.is_open()) {
219 startMultiFile_(multiFileName_);
220 }
221
222 curWriter_ = std::make_shared<VtkWriter>(gridView_, Dune::VTK::conforming);
223 curTime_ = t;
224 curOutFileName_ = fileName_();
225 }
226
233 template <class Scalar, int nComp>
234 Dune::BlockVector<Dune::FieldVector<Scalar, nComp> > *allocateManagedBuffer(int nEntities)
235 {
236 using VectorField = Dune::BlockVector<Dune::FieldVector<Scalar, nComp> >;
237
238 std::shared_ptr<ManagedVectorField_<VectorField> > vfs =
239 std::make_shared<ManagedVectorField_<VectorField> >(nEntities);
240 managedObjects_.push_back(vfs);
241 return &(vfs->vf);
242 }
243
244 template <class Scalar>
245 Dune::BlockVector<Dune::FieldVector<Scalar, 1> > *allocateManagedBuffer(int nEntities)
246 { return allocateManagedBuffer<Scalar, 1>(nEntities); }
247 Dune::BlockVector<Dune::FieldVector<double, 1> > *allocateManagedBuffer(int nEntities)
248 { return allocateManagedBuffer<double, 1>(nEntities); }
249
267 template <class DataBuffer>
268 void attachVertexData(DataBuffer &buf, std::string name, int nComps=1)
269 {
270 sanitizeBuffer_(buf, nComps);
271
272 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
274 FunctionPtr fnPtr(new VtkFn(name,
275 gridView_,
276 vertexMapper_,
277 buf,
278 /*codim=*/dim,
279 nComps));
280 curWriter_->addVertexData(fnPtr);
281 }
282
298 template <class DataBuffer>
299 void attachCellData(DataBuffer &buf, std::string name, int nComps = 1)
300 {
301 sanitizeBuffer_(buf, nComps);
302
303 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
305 FunctionPtr fnPtr(new VtkFn(name,
306 gridView_,
307 elementMapper_,
308 buf,
309 /*codim=*/0,
310 nComps));
311 curWriter_->addCellData(fnPtr);
312 }
313
321 template <class DataBuffer>
322 void attachDofData(DataBuffer &buf, std::string name, bool isVertexData, int nComps = 1)
323 {
324 if (isVertexData)
325 attachVertexData(buf, name, nComps);
326 else
327 attachCellData(buf, name, nComps);
328 }
329
337 void endWrite(bool onlyDiscard = false)
338 {
339 if (!onlyDiscard) {
340 ++curWriterNum_;
341 curWriter_->write(curOutFileName_.c_str(), OutputValue);
342
343 // determine name to write into the multi-file for the
344 // current time step
345 std::string fileName;
346 std::string suffix = fileSuffix_();
347 if (commSize_ == 1) {
348 fileName = curOutFileName_;
349 multiFile_.precision(16);
350 multiFile_ << " <DataSet timestep=\""
351 << curTime_
352 << "\" file=\""
353 << fileName << "." << suffix << "\"/>\n";
354 }
355 if (commSize_ > 1 && commRank_ == 0) {
356 // only the first process updates the multi-file
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";
364 }
365 }
366 }
367
368 // discard managed objects
369 while (managedObjects_.begin() != managedObjects_.end()) {
370 managedObjects_.pop_front();
371 }
372
373 // temporarily write the closing XML mumbo-jumbo to the mashup
374 // file so that the data set can be loaded even if the
375 // simulation is aborted (or not yet finished)
376 finishMultiFile_();
377 }
378
382 template <class Restarter>
383 void serialize(Restarter &res)
384 {
385 res.serializeSectionBegin("VTKMultiWriter");
386 res.serializeStream() << curWriterNum_ << "\n";
387
388 if (commRank_ == 0) {
389 size_t fileLen = 0;
390 size_t filePos = 0;
391 if (multiFile_.is_open()) {
392 // write the meta file into the restart file
393 filePos = multiFile_.tellp();
394 multiFile_.seekp(0, std::ios::end);
395 fileLen = multiFile_.tellp();
396 multiFile_.seekp(filePos);
397 }
398
399 res.serializeStream() << fileLen << " " << filePos << "\n";
400
401 if (fileLen > 0) {
402 std::ifstream multiFileIn(multiFileName_.c_str());
403 char *tmp = new char[fileLen];
404 multiFileIn.read(tmp, fileLen);
405 res.serializeStream().write(tmp, fileLen);
406 delete[] tmp;
407 }
408 }
409
410 res.serializeSectionEnd();
411 }
412
416 template <class Restarter>
417 void deserialize(Restarter &res)
418 {
419 res.deserializeSectionBegin("VTKMultiWriter");
420 res.deserializeStream() >> curWriterNum_;
421
422 if (commRank_ == 0) {
423 std::string dummy;
424 std::getline(res.deserializeStream(), dummy);
425
426 // recreate the meta file from the restart file
427 size_t filePos, fileLen;
428 res.deserializeStream() >> fileLen >> filePos;
429 std::getline(res.deserializeStream(), dummy);
430 if (multiFile_.is_open())
431 multiFile_.close();
432
433 if (fileLen > 0) {
434 multiFile_.open(multiFileName_.c_str());
435
436 char *tmp = new char[fileLen];
437 res.deserializeStream().read(tmp, fileLen);
438 multiFile_.write(tmp, fileLen);
439 delete[] tmp;
440 }
441
442 multiFile_.seekp(filePos);
443 }
444 else {
445 std::string tmp;
446 std::getline(res.deserializeStream(), tmp);
447 }
448 res.deserializeSectionEnd();
449 }
450
451private:
452 std::string fileName_()
453 {
454 std::ostringstream oss;
455 oss << simName_ << "-" << std::setw(5) << std::setfill('0') << curWriterNum_;
456 return oss.str();
457 }
458
459 std::string fileName_(int rank)
460 {
461 if (commSize_ > 1) {
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_;
467 return oss.str();
468
469 return oss.str();
470 }
471 else {
472 std::ostringstream oss;
473 oss << simName_ << "-" << std::setw(5) << std::setfill('0') << curWriterNum_;
474 return oss.str();
475 }
476 }
477
478 std::string fileSuffix_()
479 {
480 return (GridView::dimension == 1)?"vtp":"vtu";
481 }
482
483
484 void startMultiFile_(const std::string &multiFileName)
485 {
486 // only the first process writes to the multi-file
487 if (commRank_ == 0) {
488 // generate one meta vtk-file holding the individual time steps
489 multiFile_.open(multiFileName.c_str());
490 multiFile_ <<
491 "<?xml version=\"1.0\"?>\n"
492 "<VTKFile type=\"Collection\"\n"
493 " version=\"0.1\"\n"
494 " byte_order=\"LittleEndian\"\n"
495 " compressor=\"vtkZLibDataCompressor\">\n"
496 " <Collection>\n";
497 }
498 }
499
500 void finishMultiFile_()
501 {
502 // only the first process writes to the multi-file
503 if (commRank_ == 0) {
504 // make sure that we always have a working meta file
505 std::ofstream::pos_type pos = multiFile_.tellp();
506 multiFile_ <<
507 " </Collection>\n"
508 "</VTKFile>\n";
509 multiFile_.seekp(pos);
510 multiFile_.flush();
511 }
512 }
513
514 // make sure that all values can be displayed by paraview
515 template <class DataBuffer>
516 void sanitizeBuffer_(DataBuffer &b, int nComps)
517 {
518 for (unsigned int i = 0; i < b.size(); ++i) {
519 for (int j = 0; j < nComps; ++j) {
520 // set values which are too small to 0 to avoid
521 // problems with paraview
522 if (std::abs(b[i][j])
523 < std::numeric_limits<float>::min())
524 {
525 b[i][j] = 0.0;
526 }
527 }
528 }
529 }
530
532 // Trick: when ever we attach some data we need to copy the
533 // vector field (that's because Dune::VTKWriter is not
534 // able to write fields one at a time and using
535 // VTKWriter::add*Data doesn't copy the data's
536 // representation so that once we arrive at the point
537 // where we want to write the data to disk, it might
538 // exist anymore). The problem we encounter there is
539 // that add*Data accepts arbitrary types as vector
540 // fields, but we need a single type for the linked list
541 // which keeps track of the data added. The trick we use
542 // here is to define a non-template base class with a
543 // virtual destructor for the type given to the linked
544 // list and a derived template class which actually
545 // knows the type of the vector field it must delete.
546
547 class ManagedObject_
548 {
549 public:
550 virtual ~ManagedObject_()
551 {}
552 };
553
554 template <class VF>
555 class ManagedVectorField_ : public ManagedObject_
556 {
557 public:
558 ManagedVectorField_(int size)
559 : vf(size)
560 { }
561 VF vf;
562 };
563 // end hack
565
566 const GridView gridView_;
567 ElementMapper elementMapper_;
568 VertexMapper vertexMapper_;
569
570 std::string simName_;
571 std::ofstream multiFile_;
572 std::string multiFileName_;
573
574 int commSize_; // number of processes in the communicator
575 int commRank_; // rank of the current process in the communicator
576
577 std::shared_ptr<VtkWriter> curWriter_;
578 double curTime_;
579 std::string curOutFileName_;
580 int curWriterNum_;
581
582 std::list<std::shared_ptr<ManagedObject_> > managedObjects_;
583};
584}
585
586#endif
Helpers for deprecation.
Definition: adapt.hh:29
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