AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
VTKWriter.hpp
1#pragma once
2
3#include <string>
4#include <memory>
5
6#include <dune/grid/io/file/vtk/subsamplingvtkwriter.hh>
7#include <dune/grid/io/file/vtk/vtkwriter.hh>
8
9#include <amdis/Initfile.hpp>
10#include <amdis/common/Filesystem.hpp>
11#include <amdis/common/StaticSize.hpp>
12#include <amdis/common/ValueCategory.hpp>
13#include <amdis/gridfunctions/DiscreteFunction.hpp>
14#include <amdis/io/FileWriterBase.hpp>
15#include <amdis/io/VTKSequenceWriter.hpp>
16
17namespace AMDiS
18{
19 namespace Impl
20 {
21 template <class Tag> struct VTKFieldTypeImpl;
22 template <>
23 struct VTKFieldTypeImpl<tag::scalar> {
24 static const Dune::VTK::FieldInfo::Type value = Dune::VTK::FieldInfo::Type::scalar;
25 };
26 template <>
27 struct VTKFieldTypeImpl<tag::vector> {
28 static const Dune::VTK::FieldInfo::Type value = Dune::VTK::FieldInfo::Type::vector;
29 };
30 template <>
31 struct VTKFieldTypeImpl<tag::matrix> {
32 static const Dune::VTK::FieldInfo::Type value = Dune::VTK::FieldInfo::Type::tensor;
33 };
34
35 } // end namespace Impl
36
37
39
43 template <class GV, class GF>
45 : public FileWriterBase
46 {
47 using GridView = GV;
48 using Writer = Dune::VTKWriter<GridView>;
50
51 using GridFunction = GF;
52 using Range = typename GridFunction::Range;
53
54 template <class R>
55 static constexpr Dune::VTK::FieldInfo::Type VTKFieldType() {
56 return Impl::VTKFieldTypeImpl<ValueCategory_t<R>>::value;
57 }
58
59 template <class R>
60 static constexpr std::size_t VTKFieldSize() {
61 return static_size_v<R>;
62 }
63
64 public:
66 VTKWriter(std::string const& name, GridView const& gridView, GridFunction const& gridFunction)
67 : FileWriterBase(name)
68 , gridFunction_(gridFunction)
69 {
70 int m = 0;
71 Parameters::get(name + "->animation", animation_);
72 Parameters::get(name + "->mode", m);
73
74 mode_ = (m == 0)
75 ? Dune::VTK::ascii
76 : Dune::VTK::appendedraw;
77
78 auto subSampling = Parameters::get<int>(name + "->subsampling");
79 if (subSampling) {
80 using SubsamplingWriter = Dune::SubsamplingVTKWriter<GridView>;
81 vtkWriter_ = std::make_shared<SubsamplingWriter>(gridView, Dune::refinementIntervals(subSampling.value()));
82 } else {
83 vtkWriter_ = std::make_shared<Writer>(gridView);
84 }
85
86 if (animation_)
87 vtkSeqWriter_ = std::make_shared<SeqWriter>(vtkWriter_);
88
89 vtkWriter_->addVertexData(gridFunction_,
90 Dune::VTK::FieldInfo(name_, VTKFieldType<Range>(), VTKFieldSize<Range>()));
91 }
92
94 void write(AdaptInfo& adaptInfo, bool force) override
95 {
96 filesystem::create_directories(this->dir() + "/_piecefiles");
97 if (this->doWrite(adaptInfo) || force) {
98 if (animation_)
99 vtkSeqWriter_->pwrite(adaptInfo.time(), this->filename(), this->dir(), "_piecefiles", mode_);
100 else
101 vtkWriter_->pwrite(this->filename(), this->dir(), "_piecefiles", mode_);
102 }
103 }
104
105 private:
106 GridFunction gridFunction_;
107 std::shared_ptr<Writer> vtkWriter_;
108 std::shared_ptr<SeqWriter> vtkSeqWriter_;
109
110 // write .pvd if animation=true, otherwise write only .vtu
111 bool animation_ = false;
112
113 // represents VTK::OutputType: ascii, appendedraw
114 Dune::VTK::OutputType mode_ = Dune::VTK::ascii;
115 };
116
117} // end namespace AMDiS
Holds adapt parameters and infos about the problem.
Definition: AdaptInfo.hpp:26
double const & time() const
Gets time_.
Definition: AdaptInfo.hpp:464
Base class for filewriters.
Definition: FileWriterBase.hpp:39
std::string name_
Name of the data.
Definition: FileWriterBase.hpp:71
bool doWrite(AdaptInfo &adaptInfo) const
Return whether to write the current timestep or not.
Definition: FileWriterBase.cpp:24
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
class to write pvd-files which contains a list of all collected vtk-files
Definition: VTKSequenceWriter.hpp:28
Adapter for the dune-grid VTKWriter.
Definition: VTKWriter.hpp:46
VTKWriter(std::string const &name, GridView const &gridView, GridFunction const &gridFunction)
Constructor.
Definition: VTKWriter.hpp:66
void write(AdaptInfo &adaptInfo, bool force) override
Implements FileWriterInterface::write.
Definition: VTKWriter.hpp:94