3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
vtkfunction.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_VTK_FUNCTION_HH
25#define DUMUX_IO_VTK_FUNCTION_HH
26
27#include <string>
28#include <memory>
29
30#include <dune/grid/common/mcmgmapper.hh>
31#include <dune/grid/io/file/vtk/common.hh>
32#include <dune/grid/io/file/vtk/function.hh>
33
36
37namespace Dumux::Vtk {
38
47template <typename GridView, typename Mapper, typename F>
48struct VectorP0VTKFunction : Dune::VTKFunction<GridView>
49{
50 enum { dim = GridView::dimension };
51 using ctype = typename GridView::ctype;
52 using Element = typename GridView::template Codim<0>::Entity;
53
54public:
55
57 int ncomps() const final { return nComps_; }
58
60 std::string name() const final { return name_; }
61
63 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>&) const final
64 { return accessChooser_(mycomp, mapper_.index(e), IsIndexable<decltype(field_[0])>()); }
65
66#if DUNE_VERSION_GTE(DUNE_GRID, 2, 7)
68 Dumux::Vtk::Precision precision() const final
69 { return precision_; }
70#endif
71
73 VectorP0VTKFunction(const GridView& gridView,
74 const Mapper& mapper,
75 const F& field,
76 const std::string& name,
77 int nComps,
78 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
79 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
80 {
81 if (field.size()!=(unsigned int)(mapper.size()))
82 DUNE_THROW(Dune::IOError, "VectorP0VTKFunction: size mismatch between field "
83 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
84 }
85
86private:
87
89 double accessChooser_(int mycomp, int i, std::true_type) const
90 { return vectorFieldAccess_(mycomp, i, IsIndexable<decltype(field_[0][0])>()); }
91
93 double accessChooser_(int mycomp, int i, std::false_type) const { return field_[i]; }
94
96 double vectorFieldAccess_(int mycomp, int i, std::false_type) const { return field_[i][mycomp]; }
97
99 double vectorFieldAccess_(int mycomp, int i, std::true_type) const { DUNE_THROW(Dune::InvalidStateException, "Invalid field type"); }
100
101 const F& field_;
102 const std::string name_;
103 int nComps_;
104 const Mapper& mapper_;
105 Dumux::Vtk::Precision precision_;
106};
107
116template <typename GridView, typename Mapper, typename F>
117struct VectorP1VTKFunction : Dune::VTKFunction<GridView>
118{
119 enum { dim = GridView::dimension };
120 using ctype = typename GridView::ctype;
121 using Element = typename GridView::template Codim<0>::Entity;
122
123public:
124
126 int ncomps() const final { return nComps_; }
127
129 std::string name() const final { return name_; }
130
132 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
133 {
134 const unsigned int dim = Element::mydimension;
135 const unsigned int nVertices = e.subEntities(dim);
136
137 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
138 for (unsigned i = 0; i < nVertices; ++i)
139 cornerValues[i] = accessChooser_(mycomp, mapper_.subIndex(e, i, dim), IsIndexable<decltype(field_[0])>());
140
141 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
142 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
143 return interpolation.global(xi);
144 }
145
146#if DUNE_VERSION_GTE(DUNE_GRID, 2, 7)
148 Dumux::Vtk::Precision precision() const final
149 { return precision_; }
150#endif
151
152
154 VectorP1VTKFunction(const GridView& gridView,
155 const Mapper& mapper,
156 const F& field,
157 const std::string& name,
158 int nComps,
159 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
160 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
161 {
162 if (field.size()!=(unsigned int)( mapper.size() ))
163 DUNE_THROW(Dune::IOError, "VectorP1VTKFunction: size mismatch between field "
164 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
165 }
166private:
167
169 double accessChooser_(int mycomp, int i, std::true_type) const
170 { return vectorFieldAccess_(mycomp, i, IsIndexable<decltype(field_[0][0])>()); }
171
173 double accessChooser_(int mycomp, int i, std::false_type) const { return field_[i]; }
174
176 double vectorFieldAccess_(int mycomp, int i, std::false_type) const { return field_[i][mycomp]; }
177
179 double vectorFieldAccess_(int mycomp, int i, std::true_type) const { DUNE_THROW(Dune::InvalidStateException, "Invalid field type"); }
180
181 const F& field_;
182 const std::string name_;
183 int nComps_;
184 const Mapper& mapper_;
185 Dumux::Vtk::Precision precision_;
186};
187
200template <typename GridView, typename Mapper, typename F>
201struct VectorP1NonConformingVTKFunction : Dune::VTKFunction<GridView>
202{
203 enum { dim = GridView::dimension };
204 using ctype = typename GridView::ctype;
205 using Element = typename GridView::template Codim<0>::Entity;
206
207public:
208
210 int ncomps() const final { return nComps_; }
211
213 std::string name() const final { return name_; }
214
216 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
217 {
218 const unsigned int dim = Element::mydimension;
219 const unsigned int nVertices = e.subEntities(dim);
220
221 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
222 for (unsigned i = 0; i < nVertices; ++i)
223 cornerValues[i] = accessChooser_(mycomp, mapper_.index(e), i, IsIndexable<decltype(field_[0])>());
224
225 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
226 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
227 return interpolation.global(xi);
228 }
229
230#if DUNE_VERSION_GTE(DUNE_GRID, 2, 7)
232 Dumux::Vtk::Precision precision() const final
233 { return precision_; }
234#endif
235
237 VectorP1NonConformingVTKFunction(const GridView& gridView,
238 const Mapper& mapper,
239 const F& field,
240 const std::string& name,
241 int nComps,
242 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
243 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
244 {
245 if (field.size()!=(unsigned int)(mapper.size()))
246 DUNE_THROW(Dune::IOError, "VectorP1NonConformingVTKFunction: size mismatch between field "
247 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
248 }
249private:
250
252 double accessChooser_(int mycomp, int eIdx, int cornerIdx, std::true_type) const
253 { return fieldAccess_(mycomp, eIdx, cornerIdx, IsIndexable<decltype(field_[0][0])>()); }
254
256 double accessChooser_(int mycomp, int eIdx, int cornerIdx, std::false_type) const
257 { DUNE_THROW(Dune::InvalidStateException, "Invalid field type"); }
258
260 double fieldAccess_(int mycomp, int eIdx, int cornerIdx, std::false_type) const
261 { return field_[eIdx][cornerIdx]; }
262
264 double fieldAccess_(int mycomp, int eIdx, int cornerIdx, std::true_type) const
265 { return field_[eIdx][cornerIdx][mycomp]; }
266
267 const F& field_;
268 const std::string name_;
269 int nComps_;
270 const Mapper& mapper_;
271 Dumux::Vtk::Precision precision_;
272
273};
274
281template<class GridView>
282class Field
283{
284 enum { dim = GridView::dimension };
285 using ctype = typename GridView::ctype;
286 using Element = typename GridView::template Codim<0>::Entity;
287
288public:
289 // template constructor selects the right VTKFunction implementation
290 template <typename F, class Mapper>
291 Field(const GridView& gridView, const Mapper& mapper, F const& f,
292 const std::string& name, int numComp = 1, int codim = 0,
293 Dune::VTK::DataMode dm = Dune::VTK::conforming,
294 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
295 : codim_(codim)
296 {
297 if (codim == GridView::dimension)
298 {
299 if (dm == Dune::VTK::conforming)
300 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
301 else
302 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
303 }
304 else if (codim == 0)
305 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
306 else
307 DUNE_THROW(Dune::NotImplemented, "Only element or vertex quantities allowed.");
308 }
309
311 std::string name () const { return field_->name(); }
312
314 int ncomps() const { return field_->ncomps(); }
315
317 Dumux::Vtk::Precision precision() const
318 {
319#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
320 return Dumux::Vtk::Precision::float32;
321#else
322 return field_->precision();
323#endif
324 }
325
327 int codim() const { return codim_; }
328
330 double evaluate(int mycomp,
331 const Element &element,
332 const Dune::FieldVector< ctype, dim > &xi) const
333 { return field_->evaluate(mycomp, element, xi); }
334
336 std::shared_ptr<const Dune::VTKFunction<GridView>> get() const
337 { return field_; }
338
339private:
340 int codim_;
341 // can point to anything fulfilling the VTKFunction interface
342 std::shared_ptr<Dune::VTKFunction<GridView>> field_;
343};
344
345} // end namespace Dumux::Vtk
346
347#endif
Vtk output precision options available in Dumux.
typename Dune::is_indexable< T, I > IsIndexable
We define our own is_indexable type in order to avoid several version checks throughout dumux....
Definition: typetraits.hh:48
Definition: vtkfunction.hh:37
a VTK function that supports both scalar and vector values for each element
Definition: vtkfunction.hh:49
typename GridView::ctype ctype
Definition: vtkfunction.hh:51
typename GridView::template Codim< 0 >::Entity Element
Definition: vtkfunction.hh:52
std::string name() const final
get name
Definition: vtkfunction.hh:60
@ dim
Definition: vtkfunction.hh:50
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: vtkfunction.hh:73
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &) const final
evaluate
Definition: vtkfunction.hh:63
int ncomps() const final
return number of components
Definition: vtkfunction.hh:57
a VTK function that supports both scalar and vector values for each vertex
Definition: vtkfunction.hh:118
@ dim
Definition: vtkfunction.hh:119
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: vtkfunction.hh:132
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: vtkfunction.hh:154
typename GridView::template Codim< 0 >::Entity Element
Definition: vtkfunction.hh:121
typename GridView::ctype ctype
Definition: vtkfunction.hh:120
std::string name() const final
get name
Definition: vtkfunction.hh:129
int ncomps() const final
return number of components
Definition: vtkfunction.hh:126
A VTK function that supports both scalar and vector values for each vertex. This expects the data to ...
Definition: vtkfunction.hh:202
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: vtkfunction.hh:216
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: vtkfunction.hh:237
@ dim
Definition: vtkfunction.hh:203
typename GridView::template Codim< 0 >::Entity Element
Definition: vtkfunction.hh:205
int ncomps() const final
return number of components
Definition: vtkfunction.hh:210
typename GridView::ctype ctype
Definition: vtkfunction.hh:204
std::string name() const final
get name
Definition: vtkfunction.hh:213
struct that can hold any field that fulfills the VTKFunction interface
Definition: vtkfunction.hh:283
double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
element-local evaluation of the field
Definition: vtkfunction.hh:330
std::string name() const
return the name of this field
Definition: vtkfunction.hh:311
int ncomps() const
return the number of components of this field
Definition: vtkfunction.hh:314
std::shared_ptr< const Dune::VTKFunction< GridView > > get() const
returns the underlying vtk function
Definition: vtkfunction.hh:336
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: vtkfunction.hh:291
Dumux::Vtk::Precision precision() const
return the precision of this field
Definition: vtkfunction.hh:317
int codim() const
codimension of the entities on which the field values live
Definition: vtkfunction.hh:327