24#ifndef DUMUX_IO_VTK_FUNCTION_HH
25#define DUMUX_IO_VTK_FUNCTION_HH
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>
47template <
typename Gr
idView,
typename Mapper,
typename F>
50 enum {
dim = GridView::dimension };
51 using ctype =
typename GridView::ctype;
52 using Element =
typename GridView::template Codim<0>::Entity;
57 int ncomps() const final {
return nComps_; }
60 std::string
name() const final {
return name_; }
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])>()); }
66#if DUNE_VERSION_GTE(DUNE_GRID, 2, 7)
68 Dumux::Vtk::Precision precision() const final
69 {
return precision_; }
76 const std::string&
name,
78 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
79 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(precision)
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() <<
")");
89 double accessChooser_(
int mycomp,
int i, std::true_type)
const
90 {
return vectorFieldAccess_(mycomp, i,
IsIndexable<
decltype(field_[0][0])>()); }
93 double accessChooser_(
int mycomp,
int i, std::false_type)
const {
return field_[i]; }
96 double vectorFieldAccess_(
int mycomp,
int i, std::false_type)
const {
return field_[i][mycomp]; }
99 double vectorFieldAccess_(
int mycomp,
int i, std::true_type)
const { DUNE_THROW(Dune::InvalidStateException,
"Invalid field type"); }
102 const std::string name_;
104 const Mapper& mapper_;
105 Dumux::Vtk::Precision precision_;
116template <
typename Gr
idView,
typename Mapper,
typename F>
119 enum {
dim = GridView::dimension };
120 using ctype =
typename GridView::ctype;
121 using Element =
typename GridView::template Codim<0>::Entity;
126 int ncomps() const final {
return nComps_; }
129 std::string
name() const final {
return name_; }
132 double evaluate(
int mycomp,
const Element& e,
const Dune::FieldVector<ctype, dim>& xi)
const final
134 const unsigned int dim = Element::mydimension;
135 const unsigned int nVertices = e.subEntities(
dim);
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])>());
142 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
143 return interpolation.global(xi);
146#if DUNE_VERSION_GTE(DUNE_GRID, 2, 7)
148 Dumux::Vtk::Precision precision() const final
149 {
return precision_; }
155 const Mapper& mapper,
157 const std::string&
name,
159 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
160 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(precision)
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() <<
")");
169 double accessChooser_(
int mycomp,
int i, std::true_type)
const
170 {
return vectorFieldAccess_(mycomp, i,
IsIndexable<
decltype(field_[0][0])>()); }
173 double accessChooser_(
int mycomp,
int i, std::false_type)
const {
return field_[i]; }
176 double vectorFieldAccess_(
int mycomp,
int i, std::false_type)
const {
return field_[i][mycomp]; }
179 double vectorFieldAccess_(
int mycomp,
int i, std::true_type)
const { DUNE_THROW(Dune::InvalidStateException,
"Invalid field type"); }
182 const std::string name_;
184 const Mapper& mapper_;
185 Dumux::Vtk::Precision precision_;
200template <
typename Gr
idView,
typename Mapper,
typename F>
203 enum {
dim = GridView::dimension };
204 using ctype =
typename GridView::ctype;
205 using Element =
typename GridView::template Codim<0>::Entity;
210 int ncomps() const final {
return nComps_; }
213 std::string
name() const final {
return name_; }
216 double evaluate(
int mycomp,
const Element& e,
const Dune::FieldVector<ctype, dim>& xi)
const final
218 const unsigned int dim = Element::mydimension;
219 const unsigned int nVertices = e.subEntities(
dim);
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])>());
226 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
227 return interpolation.global(xi);
230#if DUNE_VERSION_GTE(DUNE_GRID, 2, 7)
232 Dumux::Vtk::Precision precision() const final
233 {
return precision_; }
238 const Mapper& mapper,
240 const std::string&
name,
242 Dumux::Vtk::Precision precision = Dumux::Vtk::Precision::float32)
243 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(precision)
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() <<
")");
252 double accessChooser_(
int mycomp,
int eIdx,
int cornerIdx, std::true_type)
const
253 {
return fieldAccess_(mycomp, eIdx, cornerIdx,
IsIndexable<
decltype(field_[0][0])>()); }
256 double accessChooser_(
int mycomp,
int eIdx,
int cornerIdx, std::false_type)
const
257 { DUNE_THROW(Dune::InvalidStateException,
"Invalid field type"); }
260 double fieldAccess_(
int mycomp,
int eIdx,
int cornerIdx, std::false_type)
const
261 {
return field_[eIdx][cornerIdx]; }
264 double fieldAccess_(
int mycomp,
int eIdx,
int cornerIdx, std::true_type)
const
265 {
return field_[eIdx][cornerIdx][mycomp]; }
268 const std::string name_;
270 const Mapper& mapper_;
271 Dumux::Vtk::Precision precision_;
281template<
class Gr
idView>
284 enum { dim = GridView::dimension };
285 using ctype =
typename GridView::ctype;
286 using Element =
typename GridView::template Codim<0>::Entity;
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)
297 if (
codim == GridView::dimension)
299 if (dm == Dune::VTK::conforming)
300 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
302 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
305 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
307 DUNE_THROW(Dune::NotImplemented,
"Only element or vertex quantities allowed.");
311 std::string
name ()
const {
return field_->name(); }
314 int ncomps()
const {
return field_->ncomps(); }
319#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
320 return Dumux::Vtk::Precision::float32;
322 return field_->precision();
327 int codim()
const {
return codim_; }
331 const Element &element,
332 const Dune::FieldVector< ctype, dim > &xi)
const
333 {
return field_->evaluate(mycomp, element, xi); }
336 std::shared_ptr<const Dune::VTKFunction<GridView>>
get()
const
342 std::shared_ptr<Dune::VTKFunction<GridView>> field_;
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