version 0.4.1
python/io/vtkwriter.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
11#include <dune/python/common/typeregistry.hh>
12#include <dune/python/pybind11/eigen.h>
13#include <dune/python/pybind11/functional.h>
14#include <dune/python/pybind11/pybind11.h>
15#include <dune/python/pybind11/stl.h>
16#include <dune/python/pybind11/stl_bind.h>
17#include <dune/vtk/vtkwriter.hh>
18
22
23namespace Ikarus::Python {
24
49template <class Writer, class... options>
50void registerVtkWriter(pybind11::handle scope, pybind11::class_<Writer, options...> cls) {
51 using pybind11::operator""_a;
52
54
55 using Assembler = typename Writer::Assembler;
56 using FE = typename Writer::FEType;
57 using GridView = typename Writer::GridView;
58 using VirtualizedGF = Dune::Vtk::Function<GridView>;
59
60 cls.def(pybind11::init(
61 [](std::shared_ptr<Assembler> assembler, Dune::Vtk::FormatTypes format, Dune::Vtk::DataTypes datatype,
62 Dune::Vtk::DataTypes headertype) { return new Writer(assembler, format, datatype, headertype); }),
63 pybind11::arg("assembler"), pybind11::arg("format") = Dune::Vtk::FormatTypes::BINARY,
64 pybind11::arg("datatype") = Dune::Vtk::DataTypes::FLOAT32,
65 pybind11::arg("headertype") = Dune::Vtk::DataTypes::UINT32, pybind11::keep_alive<1, 2>());
66
67 cls.def("setFormat", &Writer::setFormat);
68 cls.def("setDatatype", &Writer::setDatatype);
69 cls.def("setHeadertype", &Writer::setHeadertype);
70
71 cls.def(
72 "addAllResults", [](Writer& self, DataTag dataTag = Vtk::DataTag::asPointData) { self.addAllResults(dataTag); },
73 pybind11::arg("dataTag") = Vtk::DataTag::asPointData);
74
75 auto addResultImpl = [](Writer& self, const std::string& resType, auto type) {
76 bool success = false;
77 Dune::Hybrid::forEach(typename FE::SupportedResultTypes(), [&]<typename RT>(RT i) {
78 if (resType == toString(i)) {
79 success = true;
80 self.template addResult<RT::template Rebind>(type);
81 }
82 });
83 if (not success)
84 DUNE_THROW(Dune::NotImplemented, "Element " + Dune::className<FE>() + " doesn't support ResultType " + resType);
85 };
86
87 cls.def(
88 "addResult",
89 [&](Writer& self, const std::string& resType, DataTag tag = Vtk::DataTag::asPointData) {
90 addResultImpl(self, resType, tag);
91 },
92 pybind11::arg("resType"), pybind11::arg("tag") = Vtk::DataTag::asPointData);
93}
94
95} // namespace Ikarus::Python
Definition of the LinearElastic class for finite element mechanics computations.
void registerVtkWriter(pybind11::handle scope, pybind11::class_< Writer, options... > cls)
Register Python bindings for a VtkWriter class. .
Definition: python/io/vtkwriter.hh:50
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:82
constexpr std::string toString(DBCOption _e)
Definition: dirichletbcenforcement.hh:7
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: flatassembler.hh:21
FE class is a base class for all finite elements.
Definition: febase.hh:79
decltype(std::tuple_cat(computeSupportedResultTypes< Skills...< PreFE, typename PreFE::template FE< Skills... ... > > >()...)) SupportedResultTypes
Type alias for the supported result types by the mixin.
Definition: mixin.hh:64
Ikarus VTK Writer for finite element results.