3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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/common/typetraits.hh>
31
32#include <dune/grid/common/mcmgmapper.hh>
33#include <dune/grid/io/file/vtk/common.hh>
34#include <dune/grid/io/file/vtk/function.hh>
35
37
38namespace Dumux::Vtk {
39
40namespace Detail {
41
42template<class Field>
43double accessEntry(const Field& f, [[maybe_unused]] int mycomp, [[maybe_unused]] int i)
44{
45 if constexpr (Dune::IsIndexable<std::decay_t<decltype(std::declval<Field>()[0])>>{})
46 {
47 if constexpr (Dune::IsIndexable<std::decay_t<decltype(std::declval<Field>()[0][0])>>{})
48 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
49 else
50 return f[i][mycomp];
51 }
52 else
53 return f[i];
54}
55
56} // end namespace Detail
57
66template <typename GridView, typename Mapper, typename F>
67struct VectorP0VTKFunction : Dune::VTKFunction<GridView>
68{
69 enum { dim = GridView::dimension };
70 using ctype = typename GridView::ctype;
71 using Element = typename GridView::template Codim<0>::Entity;
72
73public:
74
76 int ncomps() const final { return nComps_; }
77
79 std::string name() const final { return name_; }
80
82 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>&) const final
83 { return Detail::accessEntry(field_, mycomp, mapper_.index(e)); }
84
86 Dumux::Vtk::Precision precision() const final
87 { return precision_; }
88
90 VectorP0VTKFunction(const GridView& gridView,
91 const Mapper& mapper,
92 const F& field,
93 const std::string& name,
94 int nComps,
95 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
96 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
97 {
98 if (field.size()!=(unsigned int)(mapper.size()))
99 DUNE_THROW(Dune::IOError, "VectorP0VTKFunction: size mismatch between field "
100 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
101 }
102
103private:
104 const F& field_;
105 const std::string name_;
106 int nComps_;
107 const Mapper& mapper_;
108 Dumux::Vtk::Precision precision_;
109};
110
119template <typename GridView, typename Mapper, typename F>
120struct VectorP1VTKFunction : Dune::VTKFunction<GridView>
121{
122 enum { dim = GridView::dimension };
123 using ctype = typename GridView::ctype;
124 using Element = typename GridView::template Codim<0>::Entity;
125
126public:
127
129 int ncomps() const final { return nComps_; }
130
132 std::string name() const final { return name_; }
133
135 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
136 {
137 const unsigned int dim = Element::mydimension;
138 const unsigned int nVertices = e.subEntities(dim);
139
140 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
141 for (unsigned i = 0; i < nVertices; ++i)
142 cornerValues[i] = Detail::accessEntry(field_, mycomp, mapper_.subIndex(e, i, dim));
143
144 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
145 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
146 return interpolation.global(xi);
147 }
148
150 Dumux::Vtk::Precision precision() const final
151 { return precision_; }
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 const F& field_;
168 const std::string name_;
169 int nComps_;
170 const Mapper& mapper_;
171 Dumux::Vtk::Precision precision_;
172};
173
186template <typename GridView, typename Mapper, typename F>
187struct VectorP1NonConformingVTKFunction : Dune::VTKFunction<GridView>
188{
189 enum { dim = GridView::dimension };
190 using ctype = typename GridView::ctype;
191 using Element = typename GridView::template Codim<0>::Entity;
192
193public:
194
196 int ncomps() const final { return nComps_; }
197
199 std::string name() const final { return name_; }
200
202 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
203 {
204 const unsigned int dim = Element::mydimension;
205 const unsigned int nVertices = e.subEntities(dim);
206
207 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
208 for (unsigned i = 0; i < nVertices; ++i)
209 cornerValues[i] = accessEntry_(mycomp, mapper_.index(e), i);
210
211 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
212 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
213 return interpolation.global(xi);
214 }
215
217 Dumux::Vtk::Precision precision() const final
218 { return precision_; }
219
221 VectorP1NonConformingVTKFunction(const GridView& gridView,
222 const Mapper& mapper,
223 const F& field,
224 const std::string& name,
225 int nComps,
226 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
227 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
228 {
229 if (field.size()!=(unsigned int)(mapper.size()))
230 DUNE_THROW(Dune::IOError, "VectorP1NonConformingVTKFunction: size mismatch between field "
231 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
232 }
233private:
235 double accessEntry_([[maybe_unused]] int mycomp, [[maybe_unused]] int eIdx, [[maybe_unused]] int cornerIdx) const
236 {
237 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0])>{})
238 {
239 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0][0])>{})
240 return field_[eIdx][cornerIdx][mycomp];
241 else
242 return field_[eIdx][cornerIdx];
243 }
244 else
245 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
246 }
247
248 const F& field_;
249 const std::string name_;
250 int nComps_;
251 const Mapper& mapper_;
252 Dumux::Vtk::Precision precision_;
253
254};
255
262template<class GridView>
263class Field
264{
265 enum { dim = GridView::dimension };
266 using ctype = typename GridView::ctype;
267 using Element = typename GridView::template Codim<0>::Entity;
268
269public:
270 // template constructor selects the right VTKFunction implementation
271 template <typename F, class Mapper>
272 Field(const GridView& gridView, const Mapper& mapper, F const& f,
273 const std::string& name, int numComp = 1, int codim = 0,
274 Dune::VTK::DataMode dm = Dune::VTK::conforming,
275 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
276 : codim_(codim)
277 {
278 if (codim == GridView::dimension)
279 {
280 if (dm == Dune::VTK::conforming)
281 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
282 else
283 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
284 }
285 else if (codim == 0)
286 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
287 else
288 DUNE_THROW(Dune::NotImplemented, "Only element or vertex quantities allowed.");
289 }
290
292 std::string name () const { return field_->name(); }
293
295 int ncomps() const { return field_->ncomps(); }
296
298 Dumux::Vtk::Precision precision() const
299 { return field_->precision(); }
300
302 int codim() const { return codim_; }
303
305 double evaluate(int mycomp,
306 const Element &element,
307 const Dune::FieldVector< ctype, dim > &xi) const
308 { return field_->evaluate(mycomp, element, xi); }
309
311 std::shared_ptr<const Dune::VTKFunction<GridView>> get() const
312 { return field_; }
313
314private:
315 int codim_;
316 // can point to anything fulfilling the VTKFunction interface
317 std::shared_ptr<Dune::VTKFunction<GridView>> field_;
318};
319
320} // end namespace Dumux::Vtk
321
322#endif
Vtk output precision options available in Dumux.
Definition: fieldtype.hh:27
double accessEntry(const Field &f, int mycomp, int i)
Definition: function.hh:43
a VTK function that supports both scalar and vector values for each element
Definition: function.hh:68
typename GridView::ctype ctype
Definition: function.hh:70
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:71
std::string name() const final
get name
Definition: function.hh:79
@ dim
Definition: function.hh:69
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:90
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &) const final
evaluate
Definition: function.hh:82
int ncomps() const final
return number of components
Definition: function.hh:76
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:86
a VTK function that supports both scalar and vector values for each vertex
Definition: function.hh:121
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:150
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: function.hh:135
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:154
@ dim
Definition: function.hh:122
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:124
typename GridView::ctype ctype
Definition: function.hh:123
std::string name() const final
get name
Definition: function.hh:132
int ncomps() const final
return number of components
Definition: function.hh:129
A VTK function that supports both scalar and vector values for each vertex. This expects the data to ...
Definition: function.hh:188
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: function.hh:202
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:221
typename GridView::template Codim< 0 >::Entity Element
Definition: function.hh:191
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: function.hh:217
int ncomps() const final
return number of components
Definition: function.hh:196
typename GridView::ctype ctype
Definition: function.hh:190
std::string name() const final
get name
Definition: function.hh:199
struct that can hold any field that fulfills the VTKFunction interface
Definition: function.hh:264
double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
element-local evaluation of the field
Definition: function.hh:305
std::string name() const
return the name of this field
Definition: function.hh:292
int ncomps() const
return the number of components of this field
Definition: function.hh:295
std::shared_ptr< const Dune::VTKFunction< GridView > > get() const
returns the underlying vtk function
Definition: function.hh:311
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:272
Dumux::Vtk::Precision precision() const
return the precision of this field
Definition: function.hh:298
int codim() const
codimension of the entities on which the field values live
Definition: function.hh:302