Dune Core Modules (2.11.0)

types.hh
1#pragma once
2
3#include <cstdint>
4#include <map>
5#include <string>
6#include <vector>
7
9#include <dune/common/typelist.hh>
11#include <dune/geometry/type.hh>
13#include <dune/vtk/utility/arguments.hh>
15
16namespace Dune::Vtk
17{
19 enum class FormatTypes {
20 UNKNOWN = 0,
21 ASCII = 1<<0,
22 BINARY = 1<<1,
23 COMPRESSED = 1<<2,
24 APPENDED = BINARY | COMPRESSED
25 };
26 std::string to_string (Vtk::FormatTypes);
27 inline const auto formatTypesList = {FormatTypes::ASCII, FormatTypes::BINARY, FormatTypes::COMPRESSED, FormatTypes::APPENDED};
28
30 Vtk::FormatTypes formatTypeOf(Dune::VTK::OutputType);
31
32
34 enum class RangeTypes {
35 UNSPECIFIED, //< The output components are not restricted
36 AUTO, //< Detect the category automatically from number of components
37 SCALAR, //< Use exactly 1 component
38 VECTOR, //< Use exactly 3 components
39 TENSOR //< Use exactly 9 components
40 };
41 std::string to_string (Vtk::RangeTypes);
42 inline const auto rangeTypesList = {RangeTypes::UNSPECIFIED, RangeTypes::AUTO, RangeTypes::SCALAR, RangeTypes::VECTOR, RangeTypes::TENSOR};
43
44 // Map a dune-grid FieldInfo::Type to ValueTypes
45 Vtk::RangeTypes rangeTypeOf (Dune::VTK::FieldInfo::Type);
46
47 // Map a number of components to a corresponding value type
48 Vtk::RangeTypes rangeTypeOf (int ncomps);
49
50
51 enum class DataTypes {
52 UNKNOWN = 0,
53 INT8, UINT8,
54 INT16, UINT16,
55 INT32, UINT32,
56 INT64, UINT64,
57 FLOAT32 = 32,
58 FLOAT64 = 64
59 };
60 std::string to_string (Vtk::DataTypes);
61 inline const auto dataTypesLists = {
62 DataTypes::UNKNOWN,
63 DataTypes::INT8, DataTypes::UINT8,
64 DataTypes::INT16, DataTypes::UINT16,
65 DataTypes::INT32, DataTypes::UINT32,
66 DataTypes::INT64, DataTypes::UINT64,
67 DataTypes::FLOAT32, DataTypes::FLOAT64
68 };
69
70 // Map a dune-grid Precision type to DataTypes
71 Vtk::DataTypes dataTypeOf (Dune::VTK::Precision);
72
73 // Map a string to DataTypes
74 Vtk::DataTypes dataTypeOf (std::string);
75
76 // Map the field_type of T to DataTypes
77 template <class T>
78 Vtk::DataTypes dataTypeOf ()
79 {
80 using F = typename FieldTraits<T>::field_type;
81 if constexpr (std::is_same_v<F, std::int8_t>) { return Vtk::DataTypes::INT8; }
82 if constexpr (std::is_same_v<F, std::uint8_t>) { return Vtk::DataTypes::UINT8; }
83 if constexpr (std::is_same_v<F, std::int16_t>) { return Vtk::DataTypes::INT16; }
84 if constexpr (std::is_same_v<F, std::uint16_t>) { return Vtk::DataTypes::UINT16; }
85 if constexpr (std::is_same_v<F, std::int32_t>) { return Vtk::DataTypes::INT32; }
86 if constexpr (std::is_same_v<F, std::uint32_t>) { return Vtk::DataTypes::UINT32; }
87 if constexpr (std::is_same_v<F, std::int64_t>) { return Vtk::DataTypes::INT64; }
88 if constexpr (std::is_same_v<F, std::uint64_t>) { return Vtk::DataTypes::UINT64; }
89 if constexpr (std::is_same_v<F, float>) { return Vtk::DataTypes::FLOAT32; }
90 if constexpr (std::is_same_v<F, double>) { return Vtk::DataTypes::FLOAT64; }
91 if constexpr (std::is_same_v<F, long double>) { return Vtk::DataTypes::FLOAT64; }
92 return Vtk::DataTypes::UNKNOWN;
93 }
94
95 template <class> struct NoConstraint : std::true_type {};
96
98 template <template <class> class C = NoConstraint, class Caller>
99 void mapDataTypes (Vtk::DataTypes t, Caller caller)
100 {
101 switch (t) {
102 case DataTypes::INT8: if constexpr(C<std::int8_t>::value) caller(MetaType<std::int8_t>{}); break;
103 case DataTypes::UINT8: if constexpr(C<std::uint8_t>::value) caller(MetaType<std::uint8_t>{}); break;
104 case DataTypes::INT16: if constexpr(C<std::int16_t>::value) caller(MetaType<std::int16_t>{}); break;
105 case DataTypes::UINT16: if constexpr(C<std::uint16_t>::value) caller(MetaType<std::uint16_t>{}); break;
106 case DataTypes::INT32: if constexpr(C<std::int32_t>::value) caller(MetaType<std::int32_t>{}); break;
107 case DataTypes::UINT32: if constexpr(C<std::uint32_t>::value) caller(MetaType<std::uint32_t>{}); break;
108 case DataTypes::INT64: if constexpr(C<std::int64_t>::value) caller(MetaType<std::int64_t>{}); break;
109 case DataTypes::UINT64: if constexpr(C<std::uint64_t>::value) caller(MetaType<std::uint64_t>{}); break;
110 case DataTypes::FLOAT32: if constexpr(C<float>::value) caller(MetaType<float>{}); break;
111 case DataTypes::FLOAT64: if constexpr(C<double>::value) caller(MetaType<double>{}); break;
112 default:
113 VTK_ASSERT_MSG(false, "Unsupported type " + to_string(t));
114 break;
115 }
116 }
117
119 template <template <class> class Constraint1 = NoConstraint,
120 template <class> class Constraint2 = NoConstraint,
121 class Caller>
122 void mapDataTypes (Vtk::DataTypes t1, Vtk::DataTypes t2, Caller caller)
123 {
124 mapDataTypes<Constraint1>(t1, [&](auto type1) {
125 mapDataTypes<Constraint2>(t2, [&](auto type2) {
126 caller(type1, type2);
127 });
128 });
129 }
130
132 template <template <class> class Constraint1 = NoConstraint,
133 template <class> class Constraint2 = NoConstraint,
134 template <class> class Constraint3 = NoConstraint,
135 class Caller>
136 void mapDataTypes (Vtk::DataTypes t1, Vtk::DataTypes t2, Vtk::DataTypes t3, Caller caller)
137 {
138 mapDataTypes<Constraint1>(t1, [&](auto type1) {
139 mapDataTypes<Constraint2>(t2, [&](auto type2) {
140 mapDataTypes<Constraint3>(t3, [&](auto type3) {
141 caller(type1, type2, type3);
142 });
143 });
144 });
145 }
146
147
148 enum class CompressorTypes {
149 NONE = 0,
150 ZLIB,
151 LZ4,
152 LZMA
153 };
154 std::string to_string (CompressorTypes);
155
156
158 struct CellType
159 {
160 enum Parametrization {
161 LINEAR = 1,
162 QUADRATIC = 2,
163 MULTI_QUADRATIC = 3,
164 LAGRANGE = 4,
165 };
166
167 enum Type : std::uint8_t {
168 UNKNOWN = 0,
169 // Linear VTK cell types
170 VERTEX = 1,
171 /* POLY_VERTEX = 2, // not supported */
172 LINE = 3,
173 /* POLY_LINE = 4, // not supported */
174 TRIANGLE = 5,
175 /* TRIANGLE_STRIP = 6, // not supported */
176 POLYGON = 7,
177 /* PIXEL = 8, // not supported */
178 QUAD = 9,
179 TETRA = 10,
180 /* VOXEL = 11, // not supported */
181 HEXAHEDRON = 12,
182 WEDGE = 13,
183 PYRAMID = 14,
184 // Quadratic VTK cell types
185 QUADRATIC_EDGE = 21,
186 QUADRATIC_TRIANGLE = 22,
187 QUADRATIC_QUAD = 23,
188 QUADRATIC_TETRA = 24,
189 QUADRATIC_HEXAHEDRON = 25,
190 QUADRATIC_WEDGE = 26,
191 QUADRATIC_PYRAMID = 27,
192 // Multi-quadratic cell types
193 BIQUADRATIC_QUAD = 28,
194 TRIQUADRATIC_HEXAHEDRON = 29,
195 BIQUADRATIC_TRIANGLE = 34,
196 TRIQUADRATIC_PYRAMID = 37,
197 // Arbitrary order Lagrange elements
198 LAGRANGE_CURVE = 68,
199 LAGRANGE_TRIANGLE = 69,
200 LAGRANGE_QUADRILATERAL= 70,
201 LAGRANGE_TETRAHEDRON = 71,
202 LAGRANGE_HEXAHEDRON = 72,
203 LAGRANGE_WEDGE = 73,
204 LAGRANGE_PYRAMID = 74,
205 };
206
207 public:
208 CellType (GeometryType const& t, Parametrization = LINEAR);
209
211 std::uint8_t type () const
212 {
213 return type_;
214 }
215
216 int layout (int dim) const
217 {
218 return layout_[dim];
219 }
220
222 int permutation (int idx) const
223 {
224 return permutation_[idx];
225 }
226
228 int size () const
229 {
230 return permutation_.size();
231 }
232
233 bool noPermutation () const
234 {
235 return noPermutation_;
236 }
237
238 private:
239 std::uint8_t type_;
240 std::vector<int> layout_; // number of nodes for each dimension
241 std::vector<int> permutation_;
242 bool noPermutation_ = true;
243 };
244 GeometryType to_geometry (std::uint8_t);
245
246
247 class FieldInfo
248 {
249 public:
250 template <class... Args>
251 explicit FieldInfo (std::string name, Args... args)
252 : name_(std::move(name))
253 , ncomps_(getArg<int,unsigned int,long,unsigned long>(args..., 1))
254 , rangeType_(getArg<Vtk::RangeTypes>(args..., Vtk::RangeTypes::AUTO))
255 , dataType_(getArg<Vtk::DataTypes>(args..., Vtk::DataTypes::FLOAT32))
256 {
257 if (rangeType_ == Vtk::RangeTypes::AUTO)
258 rangeType_ = rangeTypeOf(ncomps_);
259 }
260
261 // Construct from dune-grid FieldInfo
262 FieldInfo (Dune::VTK::FieldInfo info)
263 : FieldInfo(info.name(), info.size(), rangeTypeOf(info.type()), dataTypeOf(info.precision()))
264 {}
265
267 std::string const& name () const
268 {
269 return name_;
270 }
271
273 int size () const
274 {
275 return ncomps_;
276 }
277
279 Vtk::RangeTypes rangeType () const
280 {
281 return rangeType_;
282 }
283
285 Vtk::DataTypes dataType () const
286 {
287 return dataType_;
288 }
289
290 private:
291 std::string name_;
292 int ncomps_;
293 Vtk::RangeTypes rangeType_;
294 Vtk::DataTypes dataType_;
295 };
296
297} // end namespace Dune::Vtk
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Descriptor struct for VTK fields.
Definition: common.hh:328
Type
VTK data type.
Definition: common.hh:333
Various macros to work with Dune module version numbers.
Common stuff for the VTKWriter.
Precision
which precision to use when writing out data to vtk files
Definition: common.hh:271
OutputType
How the bulk data should be stored in the file.
Definition: common.hh:43
Macro for wrapping error checks and throwing exceptions.
#define VTK_ASSERT_MSG(cond, text)
check if condition cond holds; otherwise, throw a VtkError with a message.
Definition: errors.hh:19
Type traits to determine the type of reals (when working with complex numbers)
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
STL namespace.
Mapping of Dune geometry types to VTK cell types.
Definition: types.hh:159
int size() const
Number of nodes.
Definition: types.hh:228
std::uint8_t type() const
Return VTK Cell type.
Definition: types.hh:211
int permutation(int idx) const
Return a permutation of Dune element vertices to conform to VTK element numbering.
Definition: types.hh:222
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Feb 14, 23:39, 2026)