6 #include <dune/grid/io/file/vtk/subsamplingvtkwriter.hh> 7 #include <dune/grid/io/file/vtk/vtkwriter.hh> 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> 21 template <
class Tag>
struct VTKFieldTypeImpl;
23 struct VTKFieldTypeImpl<tag::scalar> {
24 static const Dune::VTK::FieldInfo::Type value = Dune::VTK::FieldInfo::Type::scalar;
27 struct VTKFieldTypeImpl<tag::vector> {
28 static const Dune::VTK::FieldInfo::Type value = Dune::VTK::FieldInfo::Type::vector;
31 struct VTKFieldTypeImpl<tag::matrix> {
32 static const Dune::VTK::FieldInfo::Type value = Dune::VTK::FieldInfo::Type::tensor;
43 template <
class GV,
class GF>
48 using Writer = Dune::VTKWriter<GridView>;
51 using GridFunction = GF;
52 using Range =
typename GridFunction::Range;
55 static constexpr Dune::VTK::FieldInfo::Type VTKFieldType() {
56 return Impl::VTKFieldTypeImpl<ValueCategory_t<R>>::value;
60 static constexpr std::size_t VTKFieldSize() {
61 return static_size_v<R>;
66 VTKWriter(std::string
const& name, GridView
const& gridView, GridFunction
const& gridFunction)
68 , gridFunction_(gridFunction)
76 : Dune::VTK::appendedraw;
78 auto subSampling = Parameters::get<int>(name +
"->subsampling");
80 using SubsamplingWriter = Dune::SubsamplingVTKWriter<GridView>;
81 vtkWriter_ = std::make_shared<SubsamplingWriter>(gridView, Dune::refinementIntervals(subSampling.value()));
83 vtkWriter_ = std::make_shared<Writer>(gridView);
87 vtkSeqWriter_ = std::make_shared<SeqWriter>(vtkWriter_);
89 vtkWriter_->addVertexData(gridFunction_,
90 Dune::VTK::FieldInfo(name_, VTKFieldType<Range>(), VTKFieldSize<Range>()));
96 filesystem::create_directories(this->dir() +
"/_piecefiles");
97 if (force || this->doWrite(adaptInfo)) {
99 vtkSeqWriter_->pwrite(adaptInfo.
time(), this->filename(), this->dir(),
"_piecefiles", mode_);
101 vtkWriter_->pwrite(this->filename(), this->dir(),
"_piecefiles", mode_);
106 GridFunction gridFunction_;
107 std::shared_ptr<Writer> vtkWriter_;
108 std::shared_ptr<SeqWriter> vtkSeqWriter_;
111 bool animation_ =
false;
114 Dune::VTK::OutputType mode_ = Dune::VTK::ascii;
Adapter for the dune-grid VTKWriter.
Definition: VTKWriter.hpp:44
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
VTKWriter(std::string const &name, GridView const &gridView, GridFunction const &gridFunction)
Constructor.
Definition: VTKWriter.hpp:66
Base class for filewriters.
Definition: FileWriterBase.hpp:37
class to write pvd-files which contains a list of all collected vtk-files
Definition: VTKSequenceWriter.hpp:27
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
void write(AdaptInfo &adaptInfo, bool force) override
Implements FileWriterInterface::write.
Definition: VTKWriter.hpp:94
Holds adapt parameters and infos about the problem.
Definition: AdaptInfo.hpp:25
double const & time() const
Gets time_.
Definition: AdaptInfo.hpp:464