version 3.10-dev
function.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_VTK_FUNCTION_HH
13#define DUMUX_IO_VTK_FUNCTION_HH
14
15#include <string>
16#include <memory>
17
18#include <dune/common/typetraits.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/reservedvector.hh>
21
22#include <dune/grid/common/mcmgmapper.hh>
23#include <dune/grid/io/file/vtk/common.hh>
24#include <dune/grid/io/file/vtk/function.hh>
25
28
29namespace Dumux::Vtk {
30
31namespace Detail {
32
33template<class Field>
34double accessEntry(const Field& f, [[maybe_unused]] int mycomp, [[maybe_unused]] int i)
35{
36 if constexpr (Dune::IsIndexable<std::decay_t<decltype(std::declval<Field>()[0])>>{})
37 {
38 if constexpr (Dune::IsIndexable<std::decay_t<decltype(std::declval<Field>()[0][0])>>{})
39 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
40 else
41 return f[i][mycomp];
42 }
43 else
44 return f[i];
45}
46
47} // end namespace Detail
48
57template <typename GridView, typename Mapper, typename F>
58struct VectorP0VTKFunction : Dune::VTKFunction<GridView>
59{
60 enum { dim = GridView::dimension };
61 using ctype = typename GridView::ctype;
62 using Element = typename GridView::template Codim<0>::Entity;
63
64public:
65
67 int ncomps() const final { return nComps_; }
68
70 std::string name() const final { return name_; }
71
73 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>&) const final
74 { return Detail::accessEntry(field_, mycomp, mapper_.index(e)); }
75
77 Dumux::Vtk::Precision precision() const final
78 { return precision_; }
79
81 VectorP0VTKFunction(const GridView& gridView,
82 const Mapper& mapper,
83 const F& field,
84 const std::string& name,
85 int nComps,
86 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
87 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
88 {
89 if (field.size()!=(unsigned int)(mapper.size()))
90 DUNE_THROW(Dune::IOError, "VectorP0VTKFunction: size mismatch between field "
91 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
92 }
93
94private:
95 const F& field_;
96 const std::string name_;
97 int nComps_;
98 const Mapper& mapper_;
99 Dumux::Vtk::Precision precision_;
100};
101
110template <typename GridView, typename Mapper, typename F>
111struct VectorP1VTKFunction : Dune::VTKFunction<GridView>
112{
113 enum { dim = GridView::dimension };
114 using ctype = typename GridView::ctype;
115 using Element = typename GridView::template Codim<0>::Entity;
116
117public:
118
120 int ncomps() const final { return nComps_; }
121
123 std::string name() const final { return name_; }
124
126 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
127 {
128 const unsigned int dim = Element::mydimension;
129 const unsigned int nVertices = e.subEntities(dim);
130
131 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
132 for (unsigned i = 0; i < nVertices; ++i)
133 cornerValues[i] = Detail::accessEntry(field_, mycomp, mapper_.subIndex(e, i, dim));
134
135 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
136 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
137 return interpolation.global(xi);
138 }
139
141 Dumux::Vtk::Precision precision() const final
142 { return precision_; }
143
145 VectorP1VTKFunction(const GridView& gridView,
146 const Mapper& mapper,
147 const F& field,
148 const std::string& name,
149 int nComps,
150 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
151 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
152 {
153 if (field.size()!=(unsigned int)( mapper.size() ))
154 DUNE_THROW(Dune::IOError, "VectorP1VTKFunction: size mismatch between field "
155 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
156 }
157private:
158 const F& field_;
159 const std::string name_;
160 int nComps_;
161 const Mapper& mapper_;
162 Dumux::Vtk::Precision precision_;
163};
164
177template <typename GridView, typename Mapper, typename F>
178struct VectorP1NonConformingVTKFunction : Dune::VTKFunction<GridView>
179{
180 enum { dim = GridView::dimension };
181 using ctype = typename GridView::ctype;
182 using Element = typename GridView::template Codim<0>::Entity;
183
184public:
185
187 int ncomps() const final { return nComps_; }
188
190 std::string name() const final { return name_; }
191
193 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
194 {
195 const unsigned int dim = Element::mydimension;
196 const unsigned int nVertices = e.subEntities(dim);
197
198 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
199 for (unsigned i = 0; i < nVertices; ++i)
200 cornerValues[i] = accessEntry_(mycomp, mapper_.index(e), i);
201
202 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
203 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
204 return interpolation.global(xi);
205 }
206
208 Dumux::Vtk::Precision precision() const final
209 { return precision_; }
210
212 VectorP1NonConformingVTKFunction(const GridView& gridView,
213 const Mapper& mapper,
214 const F& field,
215 const std::string& name,
216 int nComps,
217 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
218 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
219 {
220 if (field.size()!=(unsigned int)(mapper.size()))
221 DUNE_THROW(Dune::IOError, "VectorP1NonConformingVTKFunction: size mismatch between field "
222 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
223 }
224private:
226 double accessEntry_([[maybe_unused]] int mycomp, [[maybe_unused]] int eIdx, [[maybe_unused]] int cornerIdx) const
227 {
228 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0])>{})
229 {
230 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0][0])>{})
231 return field_[eIdx][cornerIdx][mycomp];
232 else
233 return field_[eIdx][cornerIdx];
234 }
235 else
236 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
237 }
238
239 const F& field_;
240 const std::string name_;
241 int nComps_;
242 const Mapper& mapper_;
243 Dumux::Vtk::Precision precision_;
244
245};
246
258template <typename GridView, typename Mapper, typename F>
259struct VectorP1FaceNonConformingVTKFunction : Dune::VTKFunction<GridView>
260{
261 static constexpr int dim = GridView::dimension;
262 using ctype = typename GridView::ctype;
263 using Element = typename GridView::template Codim<0>::Entity;
264
265public:
266
268 int ncomps() const final { return nComps_; }
269
271 std::string name() const final { return name_; }
272
274 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& localPos) const final
275 {
276 const unsigned int dim = Element::mydimension;
277 const unsigned int nFaces = e.subEntities(1);
278
279 // assemble coefficient vector
280 using Coefficients = Dune::ReservedVector<ctype, 2*dim>;
281 Coefficients faceValues;
282 faceValues.resize(nFaces);
283 for (int i = 0; i < nFaces; ++i)
284 faceValues[i] = accessEntry_(mycomp, mapper_.index(e), i);
285
286 using ShapeValues = std::vector<Dune::FieldVector<ctype, 1>>;
287 ShapeValues N;
288 feCache_.get(e.type()).localBasis().evaluateFunction(localPos, N);
289 double result = 0.0;
290 for (int i = 0; i < N.size(); ++i)
291 result += faceValues[i]*N[i];
292 return result;
293 }
294
296 Dumux::Vtk::Precision precision() const final
297 { return precision_; }
298
301 const GridView& gridView,
302 const Mapper& mapper,
303 const F& field,
304 const std::string& name,
305 int nComps,
306 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32
307 )
308 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
309 {
310 if (field.size()!=(unsigned int)(mapper.size()))
311 DUNE_THROW(Dune::IOError, "VectorP1FaceNonConformingVTKFunction: size mismatch between field "
312 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
313 }
314private:
316 double accessEntry_([[maybe_unused]] int mycomp, [[maybe_unused]] int eIdx, [[maybe_unused]] int faceIdx) const
317 {
318 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0])>{})
319 {
320 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0][0])>{})
321 return field_[eIdx][faceIdx][mycomp];
322 else
323 return field_[eIdx][faceIdx];
324 }
325 else
326 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
327 }
328
329 const F& field_;
330 const std::string name_;
331 int nComps_;
332 const Mapper& mapper_;
333 Dumux::Vtk::Precision precision_;
334 NonconformingFECache<ctype, ctype, dim> feCache_ = {};
335};
336
343template<class GridView>
344class Field
345{
346 static constexpr int dim = GridView::dimension;
347 using ctype = typename GridView::ctype;
348 using Element = typename GridView::template Codim<0>::Entity;
349
350public:
351 // template constructor selects the right VTKFunction implementation
352 template <typename F, class Mapper>
353 Field(const GridView& gridView, const Mapper& mapper, F const& f,
354 const std::string& name, int numComp = 1, int codim = 0,
355 Dune::VTK::DataMode dm = Dune::VTK::conforming,
356 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
357 : codim_(codim)
358 {
359 if constexpr (dim > 1)
360 {
361 if (codim == dim)
362 {
363 if (dm == Dune::VTK::conforming)
364 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
365 else
366 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
367 }
368 else if (codim == 1)
369 {
370 if (dm == Dune::VTK::conforming)
371 DUNE_THROW(Dune::NotImplemented, "Only non-conforming output is implemented");
372 else
373 field_ = std::make_shared< VectorP1FaceNonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
374 }
375 else if (codim == 0)
376 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
377 else
378 DUNE_THROW(Dune::NotImplemented, "Only element, face or vertex quantities allowed.");
379 }
380 else
381 {
382 if (codim == dim)
383 {
384 if (dm == Dune::VTK::conforming)
385 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
386 else
387 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
388 }
389 else if (codim == 0)
390 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
391 else
392 DUNE_THROW(Dune::NotImplemented, "Only element or vertex quantities allowed.");
393 }
394 }
395
397 std::string name () const { return field_->name(); }
398
400 int ncomps() const { return field_->ncomps(); }
401
403 Dumux::Vtk::Precision precision() const
404 { return field_->precision(); }
405
407 int codim() const { return codim_; }
408
410 double evaluate(int mycomp,
411 const Element &element,
412 const Dune::FieldVector< ctype, dim > &xi) const
413 { return field_->evaluate(mycomp, element, xi); }
414
416 std::shared_ptr<const Dune::VTKFunction<GridView>> get() const
417 { return field_; }
418
419private:
420 int codim_;
421 // can point to anything fulfilling the VTKFunction interface
422 std::shared_ptr<Dune::VTKFunction<GridView>> field_;
423};
424
425} // end namespace Dumux::Vtk
426
427#endif
const FiniteElementType & get(const Dune::GeometryType &gt) const
Get local finite element for given GeometryType.
Definition: nonconformingfecache.hh:46
struct that can hold any field that fulfills the VTKFunction interface
Definition: function.hh:345
double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
element-local evaluation of the field
Definition: function.hh:410
std::string name() const
return the name of this field
Definition: function.hh:397
int ncomps() const
return the number of components of this field
Definition: function.hh:400
std::shared_ptr< const Dune::VTKFunction< GridView > > get() const
returns the underlying vtk function
Definition: function.hh:416
Field(const GridView &gridView, const Mapper &mapper, F const &f, const std::string &name, int numComp=1, int codim=0, Dune::VTK::DataMode dm=Dune::VTK::conforming, Dumux::Vtk::Precision precision=Dumux::Vtk::Precision::float32)
Definition: function.hh:353
Dumux::Vtk::Precision precision() const
return the precision of this field
Definition: function.hh:403
int codim() const
codimension of the entities on which the field values live
Definition: function.hh:407
double accessEntry(const Field &f, int mycomp, int i)
Definition: function.hh:34
Definition: fieldtype.hh:15
A finite element cache for the non-conforming FE spaces RT and CR.
Vtk output precision options available in Dumux.
a VTK function that supports both scalar and vector values for each element
Definition: function.hh:59
typename GridView::ctype ctype
Definition: function.hh:61
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:62
std::string name() const final
get name
Definition: function.hh:70
@ dim
Definition: function.hh:60
VectorP0VTKFunction(const GridView &gridView, const Mapper &mapper, const F &field, const std::string &name, int nComps, Dumux::Vtk::Precision precision=Dumux::Vtk::Precision::float32)
Constructor.
Definition: function.hh:81
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &) const final
evaluate
Definition: function.hh:73
int ncomps() const final
return number of components
Definition: function.hh:67
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:77
A VTK function that supports both scalar and vector values for each face. This expects the data to be...
Definition: function.hh:260
std::string name() const final
get name
Definition: function.hh:271
int ncomps() const final
return number of components
Definition: function.hh:268
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &localPos) const final
evaluate
Definition: function.hh:274
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:296
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:263
VectorP1FaceNonConformingVTKFunction(const GridView &gridView, const Mapper &mapper, const F &field, const std::string &name, int nComps, Dumux::Vtk::Precision precision=Dumux::Vtk::Precision::float32)
Constructor.
Definition: function.hh:300
static constexpr int dim
Definition: function.hh:261
typename GridView::ctype ctype
Definition: function.hh:262
A VTK function that supports both scalar and vector values for each vertex. This expects the data to ...
Definition: function.hh:179
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: function.hh:193
VectorP1NonConformingVTKFunction(const GridView &gridView, const Mapper &mapper, const F &field, const std::string &name, int nComps, Dumux::Vtk::Precision precision=Dumux::Vtk::Precision::float32)
Constructor.
Definition: function.hh:212
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:182
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:208
int ncomps() const final
return number of components
Definition: function.hh:187
typename GridView::ctype ctype
Definition: function.hh:181
std::string name() const final
get name
Definition: function.hh:190
a VTK function that supports both scalar and vector values for each vertex
Definition: function.hh:112
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:141
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: function.hh:126
VectorP1VTKFunction(const GridView &gridView, const Mapper &mapper, const F &field, const std::string &name, int nComps, Dumux::Vtk::Precision precision=Dumux::Vtk::Precision::float32)
Constructor.
Definition: function.hh:145
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:115
typename GridView::ctype ctype
Definition: function.hh:114
std::string name() const final
get name
Definition: function.hh:123
@ dim
Definition: function.hh:113
int ncomps() const final
return number of components
Definition: function.hh:120