3.2-git
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 "vtknestedfunction.hh"
36
37#include <dune/common/fvector.hh>
38#include <dune/istl/bvector.hh>
39
40#include <dune/grid/io/file/vtk/vtkwriter.hh>
41
43
44
45#if HAVE_MPI
46#include <mpi.h>
47#endif
48
49namespace Dumux {
59template<class GridView, Dune::VTK::OutputType OutputValue = Dune::VTK::ascii >
60class [[deprecated("Use VtkOutputModule instead!")]] VtkMultiWriter
61{
62 enum { dim = GridView::dimension };
63 using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
64 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
65public:
66 using VtkWriter = Dune::VTKWriter<GridView>;
67 VtkMultiWriter(const GridView &gridView,
68 const std::string &simName = "",
69 std::string multiFileName = "")
70 : gridView_(gridView)
71 , elementMapper_(gridView, Dune::mcmgElementLayout())
72 , vertexMapper_(gridView, Dune::mcmgVertexLayout())
73 {
74 simName_ = (simName.empty())?"sim":simName;
75 multiFileName_ = multiFileName;
76 if (multiFileName_.empty()) {
77 multiFileName_ = simName_;
78 multiFileName_ += ".pvd";
79 }
80
81 curWriterNum_ = 0;
82
83 commRank_ = 0;
84 commSize_ = 1;
85#if HAVE_MPI
86 MPI_Comm_rank(MPI_COMM_WORLD, &commRank_);
87 MPI_Comm_size(MPI_COMM_WORLD, &commSize_);
88#endif
89 }
90
92 {
93 finishMultiFile_();
94
95 if (commRank_ == 0)
96 multiFile_.close();
97 }
98
102 int curWriterNum() const
103 { return curWriterNum_; }
104
113 {
114 elementMapper_.update();
115 vertexMapper_.update();
116 }
117
122 void beginWrite(double t)
123 {
124 if (!multiFile_.is_open()) {
125 startMultiFile_(multiFileName_);
126 }
127
128 curWriter_ = std::make_shared<VtkWriter>(gridView_, Dune::VTK::conforming);
129 curTime_ = t;
130 curOutFileName_ = fileName_();
131 }
132
139 template <class Scalar, int nComp>
140 Dune::BlockVector<Dune::FieldVector<Scalar, nComp> > *allocateManagedBuffer(int nEntities)
141 {
142 using VectorField = Dune::BlockVector<Dune::FieldVector<Scalar, nComp> >;
143
144 std::shared_ptr<ManagedVectorField_<VectorField> > vfs =
145 std::make_shared<ManagedVectorField_<VectorField> >(nEntities);
146 managedObjects_.push_back(vfs);
147 return &(vfs->vf);
148 }
149
150 template <class Scalar>
151 Dune::BlockVector<Dune::FieldVector<Scalar, 1> > *allocateManagedBuffer(int nEntities)
152 { return allocateManagedBuffer<Scalar, 1>(nEntities); }
153 Dune::BlockVector<Dune::FieldVector<double, 1> > *allocateManagedBuffer(int nEntities)
154 { return allocateManagedBuffer<double, 1>(nEntities); }
155
173 template <class DataBuffer>
174 void attachVertexData(DataBuffer &buf, std::string name, int nComps=1)
175 {
176 sanitizeBuffer_(buf, nComps);
177
178 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
180 FunctionPtr fnPtr(new VtkFn(name,
181 gridView_,
182 vertexMapper_,
183 buf,
184 /*codim=*/dim,
185 nComps));
186 curWriter_->addVertexData(fnPtr);
187 }
188
204 template <class DataBuffer>
205 void attachCellData(DataBuffer &buf, std::string name, int nComps = 1)
206 {
207 sanitizeBuffer_(buf, nComps);
208
209 using FunctionPtr = std::shared_ptr<const typename VtkWriter::VTKFunction>;
211 FunctionPtr fnPtr(new VtkFn(name,
212 gridView_,
213 elementMapper_,
214 buf,
215 /*codim=*/0,
216 nComps));
217 curWriter_->addCellData(fnPtr);
218 }
219
227 template <class DataBuffer>
228 void attachDofData(DataBuffer &buf, std::string name, bool isVertexData, int nComps = 1)
229 {
230 if (isVertexData)
231 attachVertexData(buf, name, nComps);
232 else
233 attachCellData(buf, name, nComps);
234 }
235
243 void endWrite(bool onlyDiscard = false)
244 {
245 if (!onlyDiscard) {
246 ++curWriterNum_;
247 curWriter_->write(curOutFileName_.c_str(), OutputValue);
248
249 // determine name to write into the multi-file for the
250 // current time step
251 std::string fileName;
252 std::string suffix = fileSuffix_();
253 if (commSize_ == 1) {
254 fileName = curOutFileName_;
255 multiFile_.precision(16);
256 multiFile_ << " <DataSet timestep=\""
257 << curTime_
258 << "\" file=\""
259 << fileName << "." << suffix << "\"/>\n";
260 }
261 if (commSize_ > 1 && commRank_ == 0) {
262 // only the first process updates the multi-file
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";
270 }
271 }
272 }
273
274 // discard managed objects
275 while (managedObjects_.begin() != managedObjects_.end()) {
276 managedObjects_.pop_front();
277 }
278
279 // temporarily write the closing XML mumbo-jumbo to the mashup
280 // file so that the data set can be loaded even if the
281 // simulation is aborted (or not yet finished)
282 finishMultiFile_();
283 }
284
288 template <class Restarter>
289 void serialize(Restarter &res)
290 {
291 res.serializeSectionBegin("VTKMultiWriter");
292 res.serializeStream() << curWriterNum_ << "\n";
293
294 if (commRank_ == 0) {
295 size_t fileLen = 0;
296 size_t filePos = 0;
297 if (multiFile_.is_open()) {
298 // write the meta file into the restart file
299 filePos = multiFile_.tellp();
300 multiFile_.seekp(0, std::ios::end);
301 fileLen = multiFile_.tellp();
302 multiFile_.seekp(filePos);
303 }
304
305 res.serializeStream() << fileLen << " " << filePos << "\n";
306
307 if (fileLen > 0) {
308 std::ifstream multiFileIn(multiFileName_.c_str());
309 char *tmp = new char[fileLen];
310 multiFileIn.read(tmp, fileLen);
311 res.serializeStream().write(tmp, fileLen);
312 delete[] tmp;
313 }
314 }
315
316 res.serializeSectionEnd();
317 }
318
322 template <class Restarter>
323 void deserialize(Restarter &res)
324 {
325 res.deserializeSectionBegin("VTKMultiWriter");
326 res.deserializeStream() >> curWriterNum_;
327
328 if (commRank_ == 0) {
329 std::string dummy;
330 std::getline(res.deserializeStream(), dummy);
331
332 // recreate the meta file from the restart file
333 size_t filePos, fileLen;
334 res.deserializeStream() >> fileLen >> filePos;
335 std::getline(res.deserializeStream(), dummy);
336 if (multiFile_.is_open())
337 multiFile_.close();
338
339 if (fileLen > 0) {
340 multiFile_.open(multiFileName_.c_str());
341
342 char *tmp = new char[fileLen];
343 res.deserializeStream().read(tmp, fileLen);
344 multiFile_.write(tmp, fileLen);
345 delete[] tmp;
346 }
347
348 multiFile_.seekp(filePos);
349 }
350 else {
351 std::string tmp;
352 std::getline(res.deserializeStream(), tmp);
353 }
354 res.deserializeSectionEnd();
355 }
356
357private:
358 std::string fileName_()
359 {
360 std::ostringstream oss;
361 oss << simName_ << "-" << std::setw(5) << std::setfill('0') << curWriterNum_;
362 return oss.str();
363 }
364
365 std::string fileName_(int rank)
366 {
367 if (commSize_ > 1) {
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_;
373 return oss.str();
374
375 return oss.str();
376 }
377 else {
378 std::ostringstream oss;
379 oss << simName_ << "-" << std::setw(5) << std::setfill('0') << curWriterNum_;
380 return oss.str();
381 }
382 }
383
384 std::string fileSuffix_()
385 {
386 return (GridView::dimension == 1)?"vtp":"vtu";
387 }
388
389
390 void startMultiFile_(const std::string &multiFileName)
391 {
392 // only the first process writes to the multi-file
393 if (commRank_ == 0) {
394 // generate one meta vtk-file holding the individual time steps
395 multiFile_.open(multiFileName.c_str());
396 multiFile_ <<
397 "<?xml version=\"1.0\"?>\n"
398 "<VTKFile type=\"Collection\"\n"
399 " version=\"0.1\"\n"
400 " byte_order=\"LittleEndian\"\n"
401 " compressor=\"vtkZLibDataCompressor\">\n"
402 " <Collection>\n";
403 }
404 }
405
406 void finishMultiFile_()
407 {
408 // only the first process writes to the multi-file
409 if (commRank_ == 0) {
410 // make sure that we always have a working meta file
411 std::ofstream::pos_type pos = multiFile_.tellp();
412 multiFile_ <<
413 " </Collection>\n"
414 "</VTKFile>\n";
415 multiFile_.seekp(pos);
416 multiFile_.flush();
417 }
418 }
419
420 // make sure the field is well defined if running under valgrind
421 // and make sure that all values can be displayed by paraview
422 template <class DataBuffer>
423 void sanitizeBuffer_(DataBuffer &b, int nComps)
424 {
425 for (unsigned int i = 0; i < b.size(); ++i) {
426 for (int j = 0; j < nComps; ++j) {
427 Valgrind::CheckDefined(b[i][j]);
428
429 // set values which are too small to 0 to avoid
430 // problems with paraview
431 if (std::abs(b[i][j])
432 < std::numeric_limits<float>::min())
433 {
434 b[i][j] = 0.0;
435 }
436 }
437 }
438 }
439
441 // Trick: when ever we attach some data we need to copy the
442 // vector field (that's because Dune::VTKWriter is not
443 // able to write fields one at a time and using
444 // VTKWriter::add*Data doesn't copy the data's
445 // representation so that once we arrive at the point
446 // where we want to write the data to disk, it might
447 // exist anymore). The problem we encounter there is
448 // that add*Data accepts arbitrary types as vector
449 // fields, but we need a single type for the linked list
450 // which keeps track of the data added. The trick we use
451 // here is to define a non-template base class with a
452 // virtual destructor for the type given to the linked
453 // list and a derived template class which actually
454 // knows the type of the vector field it must delete.
455
456 class ManagedObject_
457 {
458 public:
459 virtual ~ManagedObject_()
460 {}
461 };
462
463 template <class VF>
464 class ManagedVectorField_ : public ManagedObject_
465 {
466 public:
467 ManagedVectorField_(int size)
468 : vf(size)
469 { }
470 VF vf;
471 };
472 // end hack
474
475 const GridView gridView_;
476 ElementMapper elementMapper_;
477 VertexMapper vertexMapper_;
478
479 std::string simName_;
480 std::ofstream multiFile_;
481 std::string multiFileName_;
482
483 int commSize_; // number of processes in the communicator
484 int commRank_; // rank of the current process in the communicator
485
486 std::shared_ptr<VtkWriter> curWriter_;
487 double curTime_;
488 std::string curOutFileName_;
489 int curWriterNum_;
490
491 std::list<std::shared_ptr<ManagedObject_> > managedObjects_;
492};
493}
494
495#endif
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: adapt.hh:29
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