version 0.4.1
controlvtkwriter.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10#include <string>
11
12#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
13#include <dune/grid/io/file/vtk/subsamplingvtkwriter.hh>
14
15#include <Eigen/Core>
16
19
20#pragma GCC diagnostic push
21#pragma GCC diagnostic ignored "-Wswitch-enum"
22
23namespace Ikarus {
24
32template <typename B>
34{
35 using Basis = B;
36 static constexpr int components = Basis::LocalView::Tree::degree() == 0 ? 1 : Basis::LocalView::Tree::degree();
37
38public:
48 ControlSubsamplingVertexVTKWriter(const Basis& basis, const Eigen::VectorXd& sol, int refinementLevels = 0)
49 : basis_{&basis},
50 vtkWriter_(basis.gridView(), Dune::refinementLevels(refinementLevels)),
51 solution_{&sol} {}
52
53 template <typename BC>
55 this->subscribe(bc, [&](ControlMessages message) { this->updateImpl(message); });
56 return *this;
57 }
58
68 auto setFieldInfo(std::string&& name, Dune::VTK::FieldInfo::Type type, std::size_t size,
69 Dune::VTK::Precision prec = Dune::VTK::Precision::float32) {
70 fieldInfo_ = Dune::VTK::FieldInfo(std::move(name), type, size, prec);
71 isFieldInfoSet_ = true;
72 }
73
79 auto setFileNamePrefix(std::string&& name) { prefixString_ = std::move(name); }
80
90 assert(isFieldInfoSet_ && "You need to call setFieldInfo first!");
91 switch (message) {
93 auto disp = Dune::Functions::makeDiscreteGlobalBasisFunction<Dune::FieldVector<double, components>>(*basis_,
94 *solution_);
95 vtkWriter_.addVertexData(disp, fieldInfo_);
96 vtkWriter_.write(prefixString_ + std::to_string(step_++));
97 } break;
98 default:
99 break; // default: do nothing when notified
100 }
101 }
102
103private:
104 const Basis* basis_;
105 Dune::SubsamplingVTKWriter<typename Basis::GridView> vtkWriter_;
106 const Eigen::VectorXd* solution_;
107 int step_{0};
108 Dune::VTK::FieldInfo fieldInfo_{"Default", Dune::VTK::FieldInfo::Type::scalar, 1};
109 std::string prefixString_{};
110 bool isFieldInfoSet_{false};
111};
112} // namespace Ikarus
113#pragma GCC diagnostic pop
Enums for observer messages.
Implementation of the observer design pattern with broadcasters.
ControlMessages
Enum class defining control-routine-related messages.
Definition: broadcastermessages.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:34
ControlSubsamplingVertexVTKWriter(const Basis &basis, const Eigen::VectorXd &sol, int refinementLevels=0)
Constructor for ControlSubsamplingVertexVTKWriter.
Definition: controlvtkwriter.hh:48
ControlSubsamplingVertexVTKWriter & subscribeTo(BC &bc)
Definition: controlvtkwriter.hh:54
auto setFileNamePrefix(std::string &&name)
Set the file name prefix for VTK files.
Definition: controlvtkwriter.hh:79
void updateImpl(ControlMessages message)
Implementation of the update method.
Definition: controlvtkwriter.hh:89
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:68
Definition: listener.hh:27
auto subscribe(Broadcaster &broadcaster, F &&f)
Function to subscribe to a broadcaster with a given function (either a lambda, std::function or funct...
Definition: listener.hh:43