3.3.0
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
39namespace Detail {
40
41template<class Field>
42double accessEntry(const Field& f, [[maybe_unused]] int mycomp, [[maybe_unused]] int i)
43{
45 {
47 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
48 else
49 return f[i][mycomp];
50 }
51 else
52 return f[i];
53}
54
55} // end namespace Detail
56
65template <typename GridView, typename Mapper, typename F>
66struct VectorP0VTKFunction : Dune::VTKFunction<GridView>
67{
68 enum { dim = GridView::dimension };
69 using ctype = typename GridView::ctype;
70 using Element = typename GridView::template Codim<0>::Entity;
71
72public:
73
75 int ncomps() const final { return nComps_; }
76
78 std::string name() const final { return name_; }
79
81 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>&) const final
82 { return Detail::accessEntry(field_, mycomp, mapper_.index(e)); }
83
85 Dumux::Vtk::Precision precision() const final
86 { return precision_; }
87
89 VectorP0VTKFunction(const GridView& gridView,
90 const Mapper& mapper,
91 const F& field,
92 const std::string& name,
93 int nComps,
94 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
95 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
96 {
97 if (field.size()!=(unsigned int)(mapper.size()))
98 DUNE_THROW(Dune::IOError, "VectorP0VTKFunction: size mismatch between field "
99 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
100 }
101
102private:
103 const F& field_;
104 const std::string name_;
105 int nComps_;
106 const Mapper& mapper_;
107 Dumux::Vtk::Precision precision_;
108};
109
118template <typename GridView, typename Mapper, typename F>
119struct VectorP1VTKFunction : Dune::VTKFunction<GridView>
120{
121 enum { dim = GridView::dimension };
122 using ctype = typename GridView::ctype;
123 using Element = typename GridView::template Codim<0>::Entity;
124
125public:
126
128 int ncomps() const final { return nComps_; }
129
131 std::string name() const final { return name_; }
132
134 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
135 {
136 const unsigned int dim = Element::mydimension;
137 const unsigned int nVertices = e.subEntities(dim);
138
139 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
140 for (unsigned i = 0; i < nVertices; ++i)
141 cornerValues[i] = Detail::accessEntry(field_, mycomp, mapper_.subIndex(e, i, dim));
142
143 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
144 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
145 return interpolation.global(xi);
146 }
147
149 Dumux::Vtk::Precision precision() const final
150 { return precision_; }
151
153 VectorP1VTKFunction(const GridView& gridView,
154 const Mapper& mapper,
155 const F& field,
156 const std::string& name,
157 int nComps,
158 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
159 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
160 {
161 if (field.size()!=(unsigned int)( mapper.size() ))
162 DUNE_THROW(Dune::IOError, "VectorP1VTKFunction: size mismatch between field "
163 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
164 }
165private:
166 const F& field_;
167 const std::string name_;
168 int nComps_;
169 const Mapper& mapper_;
170 Dumux::Vtk::Precision precision_;
171};
172
185template <typename GridView, typename Mapper, typename F>
186struct VectorP1NonConformingVTKFunction : Dune::VTKFunction<GridView>
187{
188 enum { dim = GridView::dimension };
189 using ctype = typename GridView::ctype;
190 using Element = typename GridView::template Codim<0>::Entity;
191
192public:
193
195 int ncomps() const final { return nComps_; }
196
198 std::string name() const final { return name_; }
199
201 double evaluate(int mycomp, const Element& e, const Dune::FieldVector<ctype, dim>& xi) const final
202 {
203 const unsigned int dim = Element::mydimension;
204 const unsigned int nVertices = e.subEntities(dim);
205
206 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
207 for (unsigned i = 0; i < nVertices; ++i)
208 cornerValues[i] = accessEntry_(mycomp, mapper_.index(e), i);
209
210 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
211 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
212 return interpolation.global(xi);
213 }
214
216 Dumux::Vtk::Precision precision() const final
217 { return precision_; }
218
220 VectorP1NonConformingVTKFunction(const GridView& gridView,
221 const Mapper& mapper,
222 const F& field,
223 const std::string& name,
224 int nComps,
225 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
226 : field_(field), name_(name), nComps_(nComps), mapper_(mapper), precision_(precision)
227 {
228 if (field.size()!=(unsigned int)(mapper.size()))
229 DUNE_THROW(Dune::IOError, "VectorP1NonConformingVTKFunction: size mismatch between field "
230 << name << " (" << field.size() << ") and mapper (" << mapper.size() << ")");
231 }
232private:
234 double accessEntry_([[maybe_unused]] int mycomp, [[maybe_unused]] int eIdx, [[maybe_unused]] int cornerIdx) const
235 {
236 if constexpr (IsIndexable<decltype(std::declval<F>()[0])>{})
237 {
238 if constexpr (IsIndexable<decltype(std::declval<F>()[0][0])>{})
239 return field_[eIdx][cornerIdx][mycomp];
240 else
241 return field_[eIdx][cornerIdx];
242 }
243 else
244 DUNE_THROW(Dune::InvalidStateException, "Invalid field type");
245 }
246
247 const F& field_;
248 const std::string name_;
249 int nComps_;
250 const Mapper& mapper_;
251 Dumux::Vtk::Precision precision_;
252
253};
254
261template<class GridView>
262class Field
263{
264 enum { dim = GridView::dimension };
265 using ctype = typename GridView::ctype;
266 using Element = typename GridView::template Codim<0>::Entity;
267
268public:
269 // template constructor selects the right VTKFunction implementation
270 template <typename F, class Mapper>
271 Field(const GridView& gridView, const Mapper& mapper, F const& f,
272 const std::string& name, int numComp = 1, int codim = 0,
273 Dune::VTK::DataMode dm = Dune::VTK::conforming,
274 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
275 : codim_(codim)
276 {
277 if (codim == GridView::dimension)
278 {
279 if (dm == Dune::VTK::conforming)
280 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
281 else
282 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
283 }
284 else if (codim == 0)
285 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f, name, numComp, precision);
286 else
287 DUNE_THROW(Dune::NotImplemented, "Only element or vertex quantities allowed.");
288 }
289
291 std::string name () const { return field_->name(); }
292
294 int ncomps() const { return field_->ncomps(); }
295
297 Dumux::Vtk::Precision precision() const
298 { return field_->precision(); }
299
301 int codim() const { return codim_; }
302
304 double evaluate(int mycomp,
305 const Element &element,
306 const Dune::FieldVector< ctype, dim > &xi) const
307 { return field_->evaluate(mycomp, element, xi); }
308
310 std::shared_ptr<const Dune::VTKFunction<GridView>> get() const
311 { return field_; }
312
313private:
314 int codim_;
315 // can point to anything fulfilling the VTKFunction interface
316 std::shared_ptr<Dune::VTKFunction<GridView>> field_;
317};
318
319} // end namespace Dumux::Vtk
320
321#endif
Vtk output precision options available in Dumux.
typename Dune::IsIndexable< T, I > IsIndexable
We define our own is_indexable type in order to avoid several version checks throughout dumux....
Definition: typetraits.hh:43
Definition: vtkfunction.hh:37
double accessEntry(const Field &f, int mycomp, int i)
Definition: vtkfunction.hh:42
a VTK function that supports both scalar and vector values for each element
Definition: vtkfunction.hh:67
typename GridView::ctype ctype
Definition: vtkfunction.hh:69
typename GridView::template Codim< 0 >::Entity Element
Definition: vtkfunction.hh:70
std::string name() const final
get name
Definition: vtkfunction.hh:78
@ dim
Definition: vtkfunction.hh:68
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:89
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &) const final
evaluate
Definition: vtkfunction.hh:81
int ncomps() const final
return number of components
Definition: vtkfunction.hh:75
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: vtkfunction.hh:85
a VTK function that supports both scalar and vector values for each vertex
Definition: vtkfunction.hh:120
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: vtkfunction.hh:149
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: vtkfunction.hh:134
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:153
@ dim
Definition: vtkfunction.hh:121
typename GridView::template Codim< 0 >::Entity Element
Definition: vtkfunction.hh:123
typename GridView::ctype ctype
Definition: vtkfunction.hh:122
std::string name() const final
get name
Definition: vtkfunction.hh:131
int ncomps() const final
return number of components
Definition: vtkfunction.hh:128
A VTK function that supports both scalar and vector values for each vertex. This expects the data to ...
Definition: vtkfunction.hh:187
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const final
evaluate
Definition: vtkfunction.hh:201
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:220
typename GridView::template Codim< 0 >::Entity Element
Definition: vtkfunction.hh:190
Dumux::Vtk::Precision precision() const final
get output precision for the field
Definition: vtkfunction.hh:216
int ncomps() const final
return number of components
Definition: vtkfunction.hh:195
typename GridView::ctype ctype
Definition: vtkfunction.hh:189
std::string name() const final
get name
Definition: vtkfunction.hh:198
@ dim
Definition: vtkfunction.hh:188
struct that can hold any field that fulfills the VTKFunction interface
Definition: vtkfunction.hh:263
double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
element-local evaluation of the field
Definition: vtkfunction.hh:304
std::string name() const
return the name of this field
Definition: vtkfunction.hh:291
int ncomps() const
return the number of components of this field
Definition: vtkfunction.hh:294
std::shared_ptr< const Dune::VTKFunction< GridView > > get() const
returns the underlying vtk function
Definition: vtkfunction.hh:310
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:271
Dumux::Vtk::Precision precision() const
return the precision of this field
Definition: vtkfunction.hh:297
int codim() const
codimension of the entities on which the field values live
Definition: vtkfunction.hh:301