6#ifndef DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
7#define DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
11#include <dune/common/exceptions.hh>
12#include <dune/common/fvector.hh>
14#include <dune/geometry/type.hh>
15#include <dune/geometry/referenceelements.hh>
16#include <dune/geometry/multilineargeometry.hh>
40 template<
class Gr
idView >
46 typedef typename GridView::template Codim< 0 >::Entity
Entity;
61 const Dune::FieldVector<ctype,dim>& xi)
const = 0;
64 virtual std::string
name ()
const = 0;
94 template<
typename GV,
typename V>
130 const Dune::FieldVector<ctype,dim>&)
const override
132 return v[mapper.
index(e)*ncomps_+mycomp_];
136 std::string
name ()
const override
174 if (v.size()!=(
unsigned int)(mapper.
size()*ncomps_))
175 DUNE_THROW(IOError,
"P0VTKFunction: size mismatch");
202 template<
typename GV,
typename V>
238 const Dune::FieldVector<ctype,dim>& xi)
const override
241 const unsigned int nVertices = e.subEntities(
dim);
243 std::vector<FieldVector<ctype,1> > cornerValues(nVertices);
244 for (
unsigned i=0; i<nVertices; ++i)
245 cornerValues[i] = v[mapper.
subIndex(e,i,myDim)*ncomps_+mycomp_];
248 const MultiLinearGeometry<ctype,dim,1> interpolation(e.type(), cornerValues);
249 return interpolation.global(xi);
253 std::string
name ()
const override
291 if (v.size()!=(
unsigned int)(mapper.
size()*ncomps_))
292 DUNE_THROW(IOError,
"P1VTKFunction: size mismatch");
Mapper for multiple codim and multiple geometry types.
Common stuff for the VTKWriter.
static constexpr int dimension
The dimension of the grid.
Definition: common/gridview.hh:148
Grid::ctype ctype
type used for coordinates in grid
Definition: common/gridview.hh:145
MCMGLayout mcmgElementLayout()
layout for elements (codim-0 entities)
Definition: mcmgmapper.hh:97
MCMGLayout mcmgVertexLayout()
layout for vertices (dim-0 entities)
Definition: mcmgmapper.hh:107
Include standard header files.
Definition: agrid.hh:60
Precision
which precision to use when writing out data to vtk files
Definition: common.hh:271
static constexpr int mydimension
Dimensionality of the reference element of the entity.
Definition: common/entity.hh:112
Implementation class for a multiple codim and multiple geometry type mapper.
Definition: mcmgmapper.hh:129
size_type size() const
Return total number of entities in the entity set managed by the mapper.
Definition: mcmgmapper.hh:204
Index subIndex(const typename GV::template Codim< 0 >::Entity &e, int i, unsigned int codim) const
Map subentity of codim 0 entity to starting index in array for dof block.
Definition: mcmgmapper.hh:185
Index index(const EntityType &e) const
Map entity to starting index in array for dof block.
Definition: mcmgmapper.hh:171
A base class for grid functions with any return type and dimension.
Definition: function.hh:42
virtual double evaluate(int comp, const Entity &e, const Dune::FieldVector< ctype, dim > &xi) const =0
evaluate single component comp in the entity e at local coordinates xi
GridView::ctype ctype
Definition: function.hh:44
GridView::template Codim< 0 >::Entity Entity
Definition: function.hh:46
static constexpr int dim
Definition: function.hh:45
virtual std::string name() const =0
get name
virtual VTK::Precision precision() const
get output precision for the field
Definition: function.hh:67
virtual int ncomps() const =0
virtual ~VTKFunction()
virtual destructor
Definition: function.hh:71
Take a vector and interpret it as cell data for the VTKWriter.
Definition: function.hh:97
std::string name() const override
get name
Definition: function.hh:136
Base::ctype ctype
Definition: function.hh:119
double evaluate(int, const Entity &e, const Dune::FieldVector< ctype, dim > &) const override
evaluate
Definition: function.hh:129
P0VTKFunction(const GV &gv, const V &v_, const std::string &s_, int ncomps=1, int mycomp=0, VTK::Precision prec=VTK::Precision::float32)
construct from a vector and a name
Definition: function.hh:165
VTK::Precision precision() const override
get output precision for the field
Definition: function.hh:142
Base::Entity Entity
Definition: function.hh:118
virtual ~P0VTKFunction()
destructor
Definition: function.hh:179
int ncomps() const override
return number of components
Definition: function.hh:123
Take a vector and interpret it as point data for the VTKWriter.
Definition: function.hh:205
std::string name() const override
get name
Definition: function.hh:253
virtual ~P1VTKFunction()
destructor
Definition: function.hh:296
Base::Entity Entity
Definition: function.hh:226
Base::ctype ctype
Definition: function.hh:227
VTK::Precision precision() const override
get output precision for the field
Definition: function.hh:259
double evaluate(int comp, const Entity &e, const Dune::FieldVector< ctype, dim > &xi) const override
evaluate
Definition: function.hh:237
P1VTKFunction(const GV &gv, const V &v_, const std::string &s_, int ncomps=1, int mycomp=0, VTK::Precision prec=VTK::Precision::float32)
construct from a vector and a name
Definition: function.hh:282
int ncomps() const override
return number of components
Definition: function.hh:231