AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
DuneVtkWriter.hpp
1#pragma once
2
3#include <string>
4#include <memory>
5
6#if HAVE_DUNE_VTK
7
8#include <dune/vtk/pvdwriter.hh>
9#include <dune/vtk/vtkwriter.hh>
10
11#include <amdis/Initfile.hpp>
12#include <amdis/common/Filesystem.hpp>
13#include <amdis/io/FileWriterBase.hpp>
14
15namespace AMDiS
16{
18
22 template <class GV, class GF>
23 class DuneVtkWriter
24 : public FileWriterBase
25 {
26 using GridView = GV;
27 using Writer = Dune::Vtk::VtkWriter<GridView>;
28 using SeqWriter = Dune::Vtk::PvdWriter<Writer>;
29 using GridFunction = GF;
30
31 public:
33 DuneVtkWriter(std::string const& name, GridView const& gridView, GridFunction const& gridFunction)
34 : FileWriterBase(name)
35 , gridFunction_(gridFunction)
36 {
37 int m = 0, p = 0;
38 Parameters::get(name + "->animation", animation_);
39 Parameters::get(name + "->mode", m);
40 Parameters::get(name + "->precision", p);
41
42 Dune::Vtk::FormatTypes mode =
43 m == 0 ? Dune::Vtk::FormatTypes::ASCII :
44 m == 1 ? Dune::Vtk::FormatTypes::BINARY :
45 Dune::Vtk::FormatTypes::COMPRESSED;
46
47 Dune::Vtk::DataTypes precision =
48 p == 0 ? Dune::Vtk::DataTypes::FLOAT32 :
49 Dune::Vtk::DataTypes::FLOAT64;
50
51 if (animation_) {
52 vtkSeqWriter_ = std::make_shared<SeqWriter>(gridView, mode, precision);
53 vtkSeqWriter_->addPointData(gridFunction_, this->name());
54 } else {
55 vtkWriter_ = std::make_shared<Writer>(gridView, mode, precision);
56 vtkWriter_->addPointData(gridFunction_, this->name());
57 }
58 }
59
61 void write(AdaptInfo& adaptInfo, bool force) override
62 {
63 std::string data_dir = this->dir() + "/_piecefiles";
64 filesystem::create_directories(data_dir);
65
66 std::string filename = filesystem::path({this->dir(), this->filename() + ".vtk"}).string();
67 if (this->doWrite(adaptInfo) || force) {
68 if (animation_)
69 vtkSeqWriter_->writeTimestep(adaptInfo.time(), filename, data_dir, true);
70 else
71 vtkWriter_->write(filename, data_dir);
72 }
73 }
74
75 private:
76 GridFunction gridFunction_;
77
78 std::shared_ptr<Writer> vtkWriter_;
79 std::shared_ptr<SeqWriter> vtkSeqWriter_;
80
82 bool animation_ = false;
83 };
84
85} // end namespace AMDiS
86
87#endif // HAVE_DUNE_VTK
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
constexpr bool GridFunction
GridFunction GF is a Type that has LocalFunction and provides some typedefs for Domain,...
Definition: GridFunction.hpp:72