12#ifndef DUMUX_IO_VTK_FUNCTION_HH
13#define DUMUX_IO_VTK_FUNCTION_HH
18#include <dune/common/typetraits.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/reservedvector.hh>
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>
36 if constexpr (Dune::IsIndexable<std::decay_t<decltype(std::declval<Field>()[0])>>{})
38 if constexpr (Dune::IsIndexable<std::decay_t<decltype(std::declval<Field>()[0][0])>>{})
39 DUNE_THROW(Dune::InvalidStateException,
"Invalid field type");
57template <
typename Gr
idView,
typename Mapper,
typename F>
60 enum {
dim = GridView::dimension };
61 using ctype =
typename GridView::ctype;
62 using Element =
typename GridView::template Codim<0>::Entity;
67 int ncomps() const final {
return nComps_; }
70 std::string
name() const final {
return name_; }
73 double evaluate(
int mycomp,
const Element& e,
const Dune::FieldVector<ctype, dim>&)
const final
78 {
return precision_; }
84 const std::string&
name,
86 Dumux::Vtk::Precision
precision = Dumux::Vtk::Precision::float32)
87 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(
precision)
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() <<
")");
96 const std::string name_;
98 const Mapper& mapper_;
99 Dumux::Vtk::Precision precision_;
110template <
typename Gr
idView,
typename Mapper,
typename F>
113 enum {
dim = GridView::dimension };
114 using ctype =
typename GridView::ctype;
115 using Element =
typename GridView::template Codim<0>::Entity;
120 int ncomps() const final {
return nComps_; }
123 std::string
name() const final {
return name_; }
126 double evaluate(
int mycomp,
const Element& e,
const Dune::FieldVector<ctype, dim>& xi)
const final
128 const unsigned int dim = Element::mydimension;
129 const unsigned int nVertices = e.subEntities(
dim);
131 std::vector<Dune::FieldVector<ctype, 1>> cornerValues(nVertices);
132 for (
unsigned i = 0; i < nVertices; ++i)
136 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
137 return interpolation.global(xi);
142 {
return precision_; }
146 const Mapper& mapper,
148 const std::string&
name,
150 Dumux::Vtk::Precision
precision = Dumux::Vtk::Precision::float32)
151 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(
precision)
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() <<
")");
159 const std::string name_;
161 const Mapper& mapper_;
162 Dumux::Vtk::Precision precision_;
177template <
typename Gr
idView,
typename Mapper,
typename F>
180 enum {
dim = GridView::dimension };
181 using ctype =
typename GridView::ctype;
182 using Element =
typename GridView::template Codim<0>::Entity;
187 int ncomps() const final {
return nComps_; }
190 std::string
name() const final {
return name_; }
193 double evaluate(
int mycomp,
const Element& e,
const Dune::FieldVector<ctype, dim>& xi)
const final
195 const unsigned int dim = Element::mydimension;
196 const unsigned int nVertices = e.subEntities(
dim);
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);
203 const Dune::MultiLinearGeometry<ctype, dim, 1> interpolation(e.type(), std::move(cornerValues));
204 return interpolation.global(xi);
209 {
return precision_; }
213 const Mapper& mapper,
215 const std::string&
name,
217 Dumux::Vtk::Precision
precision = Dumux::Vtk::Precision::float32)
218 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(
precision)
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() <<
")");
226 double accessEntry_([[maybe_unused]]
int mycomp, [[maybe_unused]]
int eIdx, [[maybe_unused]]
int cornerIdx)
const
228 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0])>{})
230 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0][0])>{})
231 return field_[eIdx][cornerIdx][mycomp];
233 return field_[eIdx][cornerIdx];
236 DUNE_THROW(Dune::InvalidStateException,
"Invalid field type");
240 const std::string name_;
242 const Mapper& mapper_;
243 Dumux::Vtk::Precision precision_;
258template <
typename Gr
idView,
typename Mapper,
typename F>
261 static constexpr int dim = GridView::dimension;
262 using ctype =
typename GridView::ctype;
263 using Element =
typename GridView::template Codim<0>::Entity;
268 int ncomps() const final {
return nComps_; }
271 std::string
name() const final {
return name_; }
274 double evaluate(
int mycomp,
const Element& e,
const Dune::FieldVector<ctype, dim>& localPos)
const final
276 const unsigned int dim = Element::mydimension;
277 const unsigned int nFaces = e.subEntities(1);
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);
286 using ShapeValues = std::vector<Dune::FieldVector<ctype, 1>>;
288 feCache_.
get(e.type()).localBasis().evaluateFunction(localPos, N);
290 for (
int i = 0; i < N.size(); ++i)
291 result += faceValues[i]*N[i];
297 {
return precision_; }
301 const GridView& gridView,
302 const Mapper& mapper,
304 const std::string&
name,
306 Dumux::Vtk::Precision
precision = Dumux::Vtk::Precision::float32
308 : field_(field), name_(
name), nComps_(nComps), mapper_(mapper), precision_(
precision)
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() <<
")");
316 double accessEntry_([[maybe_unused]]
int mycomp, [[maybe_unused]]
int eIdx, [[maybe_unused]]
int faceIdx)
const
318 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0])>{})
320 if constexpr (Dune::IsIndexable<decltype(std::declval<F>()[0][0])>{})
321 return field_[eIdx][faceIdx][mycomp];
323 return field_[eIdx][faceIdx];
326 DUNE_THROW(Dune::InvalidStateException,
"Invalid field type");
330 const std::string name_;
332 const Mapper& mapper_;
333 Dumux::Vtk::Precision precision_;
334 NonconformingFECache<ctype, ctype, dim> feCache_ = {};
343template<
class Gr
idView>
346 static constexpr int dim = GridView::dimension;
347 using ctype =
typename GridView::ctype;
348 using Element =
typename GridView::template Codim<0>::Entity;
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)
359 if constexpr (dim > 1)
363 if (dm == Dune::VTK::conforming)
364 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
366 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
370 if (dm == Dune::VTK::conforming)
371 DUNE_THROW(Dune::NotImplemented,
"Only non-conforming output is implemented");
373 field_ = std::make_shared< VectorP1FaceNonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
376 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
378 DUNE_THROW(Dune::NotImplemented,
"Only element, face or vertex quantities allowed.");
384 if (dm == Dune::VTK::conforming)
385 field_ = std::make_shared< VectorP1VTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
387 field_ = std::make_shared< VectorP1NonConformingVTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
390 field_ = std::make_shared< VectorP0VTKFunction<GridView, Mapper, F> >(gridView, mapper, f,
name, numComp,
precision);
392 DUNE_THROW(Dune::NotImplemented,
"Only element or vertex quantities allowed.");
397 std::string
name ()
const {
return field_->name(); }
400 int ncomps()
const {
return field_->ncomps(); }
404 {
return field_->precision(); }
407 int codim()
const {
return codim_; }
411 const Element &element,
412 const Dune::FieldVector< ctype, dim > &xi)
const
413 {
return field_->evaluate(mycomp,
element, xi); }
416 std::shared_ptr<const Dune::VTKFunction<GridView>>
get()
const
422 std::shared_ptr<Dune::VTKFunction<GridView>> field_;
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
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 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