version 0.4.1
controlvtkwriter.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10#include "observer.hh"
11#include "observermessages.hh"
12
13#include <string>
14
15#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
16#include <dune/grid/io/file/vtk/subsamplingvtkwriter.hh>
17
18#pragma GCC diagnostic push
19#pragma GCC diagnostic ignored "-Wswitch-enum"
20
21namespace Ikarus {
22
30template <typename B>
31class ControlSubsamplingVertexVTKWriter : public IObserver<ControlMessages>
32{
33 using Basis = B;
34 static constexpr int components = Basis::LocalView::Tree::degree() == 0 ? 1 : Basis::LocalView::Tree::degree();
35
36public:
46 ControlSubsamplingVertexVTKWriter(const Basis& basis, const Eigen::VectorXd& sol, int refinementLevels = 0)
47 : basis_{&basis},
48 vtkWriter_(basis.gridView(), Dune::refinementLevels(refinementLevels)),
49 solution_{&sol} {}
50
60 auto setFieldInfo(std::string&& name, Dune::VTK::FieldInfo::Type type, std::size_t size,
61 Dune::VTK::Precision prec = Dune::VTK::Precision::float32) {
62 fieldInfo_ = Dune::VTK::FieldInfo(std::move(name), type, size, prec);
63 isFieldInfoSet_ = true;
64 }
65
71 auto setFileNamePrefix(std::string&& name) { prefixString_ = std::move(name); }
72
81 void updateImpl(ControlMessages message) final {
82 assert(isFieldInfoSet_ && "You need to call setFieldInfo first!");
83 switch (message) {
85 auto disp = Dune::Functions::makeDiscreteGlobalBasisFunction<Dune::FieldVector<double, components>>(*basis_,
86 *solution_);
87 vtkWriter_.addVertexData(disp, fieldInfo_);
88 vtkWriter_.write(prefixString_ + std::to_string(step_++));
89 } break;
90 default:
91 break; // default: do nothing when notified
92 }
93 }
94
95private:
96 const Basis* basis_;
97 Dune::SubsamplingVTKWriter<typename Basis::GridView> vtkWriter_;
98 const Eigen::VectorXd* solution_;
99 int step_{0};
100 Dune::VTK::FieldInfo fieldInfo_{"Default", Dune::VTK::FieldInfo::Type::scalar, 1};
101 std::string prefixString_{};
102 bool isFieldInfoSet_{false};
103};
104} // namespace Ikarus
105#pragma GCC diagnostic pop
Enums for observer messages.
Implementation of the observer design pattern.
ControlMessages
Enum class defining control-routine-related messages.
Definition: observermessages.hh:17
Definition: assemblermanipulatorbuildingblocks.hh:22
Definition: utils/dirichletvalues.hh:30
def basis(gv, tree)
Definition: basis.py:10
ControlSubsamplingVertexVTKWriter class for writing VTK files with subsampling based on control messa...
Definition: controlvtkwriter.hh:32
ControlSubsamplingVertexVTKWriter(const Basis &basis, const Eigen::VectorXd &sol, int refinementLevels=0)
Constructor for ControlSubsamplingVertexVTKWriter.
Definition: controlvtkwriter.hh:46
void updateImpl(ControlMessages message) final
Implementation of the update method.
Definition: controlvtkwriter.hh:81
auto setFileNamePrefix(std::string &&name)
Set the file name prefix for VTK files.
Definition: controlvtkwriter.hh:71
auto setFieldInfo(std::string &&name, Dune::VTK::FieldInfo::Type type, std::size_t size, Dune::VTK::Precision prec=Dune::VTK::Precision::float32)
Set field information for the VTK file.
Definition: controlvtkwriter.hh:60
Generic observer interface for the Observer design pattern. See for a description of the design patt...
Definition: observer.hh:26