version 0.4.7
io/vtkwriter.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers ikarus@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
11#pragma once
12
13#include "vtkcontainertypes.hh"
14#include "vtkdatatag.hh"
15
16#include <dune/common/fvector.hh>
17#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
18#include <dune/grid/yaspgrid.hh>
19#include <dune/vtk/vtkwriter.hh>
20#include <dune/vtk/writers/unstructuredgridwriter.hh>
21
24
25namespace Ikarus::Vtk {
26
35template <typename AS, typename DC, typename Base>
36requires(Concepts::FlatAssembler<AS> && Concepts::DataCollector<DC>)
37struct Writer : public Base
38{
39public:
40 using Assembler = AS;
41 using GridView = typename Assembler::GridView;
42 using FERequirement = typename Assembler::FERequirement;
43 using FEContainer = typename Assembler::FEContainer;
44 using FEType = typename std::remove_cvref_t<FEContainer>::value_type;
45
46 using DataCollector = DC;
47 using VTKWriter = Base;
48
55 template <class... Args>
56 Writer(std::shared_ptr<AS> assembler, Args... args)
57 : Base(assembler->gridView(), std::forward<Args>(args)...),
58 assembler_(assembler) {}
59
67 template <typename DC_, class... Args>
69 Writer(std::shared_ptr<AS> assembler, DC_&& dc, Args... args)
70 : Base(std::forward<std::decay_t<DC_>>(dc), std::forward<Args>(args)...),
71 assembler_(assembler) {}
72
80 template <typename RF>
81 void addResultFunction(RF&& resultFunction, DataTag dataTag = DataTag::asPointData) {
82 if (dataTag == DataTag::asCellData or dataTag == DataTag::asCellAndPointData)
83 Base::addCellData(std::forward<RF>(resultFunction));
84 if (dataTag == DataTag::asPointData or dataTag == DataTag::asCellAndPointData)
85 Base::addPointData(std::forward<RF>(resultFunction));
86 }
87
94 template <template <typename, int, int> class RT, typename UserFunction = Ikarus::Impl::DefaultUserFunction>
96 void addResult(DataTag dataTag = DataTag::asPointData, UserFunction&& userFunction = {}) {
97 auto resFunction = makeResultFunction<RT>(assembler_, std::forward<UserFunction>(userFunction));
98 addResultFunction(std::move(resFunction), dataTag);
99 }
100
107 using ResultTuple = typename FEType::SupportedResultTypes;
108
109 Dune::Hybrid::forEach(ResultTuple(), [&]<typename RT>(RT i) { addResult<RT::template Rebind>(dataTag); });
110 }
111
125 template <typename Basis, typename R>
126 void addInterpolation(R&& vals, Basis&& basis, const std::string& name, DataTag dataTag = DataTag::asPointData) {
127 using Container = Impl::ResultContainer_t<std::remove_cvref_t<Basis>>;
128
129 auto gridFunction =
130 Dune::Functions::makeDiscreteGlobalBasisFunction<Container>(std::forward<Basis>(basis), std::forward<R>(vals));
131 auto fieldInfo = Dune::Vtk::FieldInfo(name, Impl::sizeOfContainer<Container>);
132
133 if (dataTag == DataTag::asCellData or dataTag == DataTag::asCellAndPointData)
134 Base::addCellData(gridFunction, fieldInfo);
135 if (dataTag == DataTag::asPointData or dataTag == DataTag::asCellAndPointData)
136 Base::addPointData(gridFunction, fieldInfo);
137 }
138
139private:
140 std::shared_ptr<Assembler> assembler_;
141};
142
148template <typename G>
149struct IsStructured : std::false_type
150{
151};
152
156template <int dim, typename Coordinates>
157struct IsStructured<Dune::YaspGrid<dim, Coordinates>> : std::true_type
158{
159};
160
166template <typename GV>
169{
172 std::conditional_t<isStructured, Dune::Vtk::YaspDataCollector<GV>, Dune::Vtk::ContinuousDataCollector<GV>>;
173
174 template <typename DC = DefaultDataCollector>
175 using DefaultVTKWriter = std::conditional_t<isStructured, Dune::Vtk::RectilinearGridWriter<typename DC::GridView, DC>,
176 Dune::Vtk::UnstructuredGridWriter<typename DC::GridView, DC>>;
177};
178
179// Class template argument deduction guides for VTK::Writer
180
181template <typename AS, class... Args>
183Writer(std::shared_ptr<AS>,
186
187template <typename AS, typename DC, class... Args, Dune::Vtk::IsDataCollector<std::decay_t<DC>> = true>
189Writer(std::shared_ptr<AS>, DC&&, Args...)
191 typename DefaultVTKWriterManager<typename AS::GridView>::template DefaultVTKWriter<std::decay_t<DC>>>;
192
193} // namespace Ikarus::Vtk
Implementation of container types for the vtk writer.
Ikarus Result Function for Stress and other finite element results.
Definition: vtkcontainertypes.hh:19
Writer(std::shared_ptr< AS >, Args...) -> Writer< AS, typename DefaultVTKWriterManager< typename AS::GridView >::DefaultDataCollector, typename DefaultVTKWriterManager< typename AS::GridView >::template DefaultVTKWriter<> >
DataTag
Tag enum indicating cell data or point data.
Definition: vtkdatatag.hh:13
Definition: utils/dirichletvalues.hh:35
def basis(gv, tree)
Definition: basis.py:10
Manages writing results using VTK, based on assembler and data collector.
Definition: io/vtkwriter.hh:38
void addAllResults(DataTag dataTag=DataTag::asPointData)
Adds all results for the given data tag.
Definition: io/vtkwriter.hh:106
Writer(std::shared_ptr< AS > assembler, DC_ &&dc, Args... args)
Constructor with assembler, data collector, and additional arguments.
Definition: io/vtkwriter.hh:69
typename std::remove_cvref_t< FEContainer >::value_type FEType
Definition: io/vtkwriter.hh:44
AS Assembler
Definition: io/vtkwriter.hh:40
void addResultFunction(RF &&resultFunction, DataTag dataTag=DataTag::asPointData)
Adds a result function for the given data tag.
Definition: io/vtkwriter.hh:81
DC DataCollector
Definition: io/vtkwriter.hh:46
Base VTKWriter
Definition: io/vtkwriter.hh:47
Writer(std::shared_ptr< AS > assembler, Args... args)
Constructor with assembler and additional arguments.
Definition: io/vtkwriter.hh:56
void addResult(DataTag dataTag=DataTag::asPointData, UserFunction &&userFunction={})
Adds a result for the given data tag.
Definition: io/vtkwriter.hh:96
void addInterpolation(R &&vals, Basis &&basis, const std::string &name, DataTag dataTag=DataTag::asPointData)
Adds interpolation data for the given basis and container.
Definition: io/vtkwriter.hh:126
typename Assembler::FEContainer FEContainer
Definition: io/vtkwriter.hh:43
typename Assembler::GridView GridView
Definition: io/vtkwriter.hh:41
typename Assembler::FERequirement FERequirement
Definition: io/vtkwriter.hh:42
Meta type to check whether a grid is structured, inherits from false_type.
Definition: io/vtkwriter.hh:150
Manages the default template parameter for the Vtk::Writer
Definition: io/vtkwriter.hh:169
static constexpr bool isStructured
Definition: io/vtkwriter.hh:170
std::conditional_t< isStructured, Dune::Vtk::RectilinearGridWriter< typename DC::GridView, DC >, Dune::Vtk::UnstructuredGridWriter< typename DC::GridView, DC > > DefaultVTKWriter
Definition: io/vtkwriter.hh:176
std::conditional_t< isStructured, Dune::Vtk::YaspDataCollector< GV >, Dune::Vtk::ContinuousDataCollector< GV > > DefaultDataCollector
Definition: io/vtkwriter.hh:172
A concept to check if a template type satisfies the ResultType requirements.
Definition: utils/concepts.hh:505
Concept representing the requirements for a FlatAssembler.A type T satisfies FlatAssembler if it prov...
Definition: utils/concepts.hh:517
Definition: utils/concepts.hh:589
Definition: utils/concepts.hh:598
Several concepts.