Dune Core Modules (2.7.1)

function.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
5#define DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
6
7#include <string>
8
11
12#include <dune/geometry/type.hh>
13#include <dune/geometry/referenceelements.hh>
14#include <dune/geometry/multilineargeometry.hh>
15
18
24namespace Dune
25{
28
30 //
31 // Base VTKFunction
32 //
33
38 template< class GridView >
40 {
41 public:
42 typedef typename GridView::ctype ctype;
43 enum { dim = GridView::dimension };
44 typedef typename GridView::template Codim< 0 >::Entity Entity;
45
48 virtual int ncomps () const = 0;
49
51
58 virtual double evaluate (int comp, const Entity& e,
59 const Dune::FieldVector<ctype,dim>& xi) const = 0;
60
62 virtual std::string name () const = 0;
63
65 virtual VTK::Precision precision() const
66 { return VTK::Precision::float32; }
67
69 virtual ~VTKFunction () {}
70 };
71
73 //
74 // P0VTKFunction
75 //
76
78
92 template<typename GV, typename V>
94 : public VTKFunction< GV >
95 {
97 typedef VTKFunction< GV > Base;
100
102 const V& v;
104 std::string s;
106 int ncomps_;
109 int mycomp_;
111 VTK::Precision prec_;
113 Mapper mapper;
114
115 public:
116 typedef typename Base::Entity Entity;
117 typedef typename Base::ctype ctype;
118 using Base::dim;
119
121 int ncomps () const override
122 {
123 return 1;
124 }
125
127 double evaluate (int, const Entity& e,
128 const Dune::FieldVector<ctype,dim>&) const override
129 {
130 return v[mapper.index(e)*ncomps_+mycomp_];
131 }
132
134 std::string name () const override
135 {
136 return s;
137 }
138
140 VTK::Precision precision() const override
141 {
142 return prec_;
143 }
144
146
163 P0VTKFunction(const GV &gv, const V &v_, const std::string &s_,
164 int ncomps=1, int mycomp=0, VTK::Precision prec = VTK::Precision::float32)
165 : v( v_ ),
166 s( s_ ),
167 ncomps_(ncomps),
168 mycomp_(mycomp),
169 prec_(prec),
170 mapper( gv, mcmgElementLayout() )
171 {
172 if (v.size()!=(unsigned int)(mapper.size()*ncomps_))
173 DUNE_THROW(IOError, "P0VTKFunction: size mismatch");
174 }
175
177 virtual ~P0VTKFunction() {}
178 };
179
181 //
182 // P1VTKFunction
183 //
184
186
200 template<typename GV, typename V>
202 : public VTKFunction< GV >
203 {
205 typedef VTKFunction< GV > Base;
208
210 const V& v;
212 std::string s;
214 int ncomps_;
217 int mycomp_;
219 VTK::Precision prec_;
221 Mapper mapper;
222
223 public:
224 typedef typename Base::Entity Entity;
225 typedef typename Base::ctype ctype;
226 using Base::dim;
227
229 int ncomps () const override
230 {
231 return 1;
232 }
233
235 double evaluate (int comp, const Entity& e,
236 const Dune::FieldVector<ctype,dim>& xi) const override
237 {
238 const unsigned int dim = Entity::mydimension;
239 const unsigned int nVertices = e.subEntities(dim);
240
241 std::vector<FieldVector<ctype,1> > cornerValues(nVertices);
242 for (unsigned i=0; i<nVertices; ++i)
243 cornerValues[i] = v[mapper.subIndex(e,i,dim)*ncomps_+mycomp_];
244
245 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
246 const MultiLinearGeometry<ctype,dim,1> interpolation(e.type(), cornerValues);
247 return interpolation.global(xi);
248 }
249
251 std::string name () const override
252 {
253 return s;
254 }
255
257 VTK::Precision precision() const override
258 {
259 return prec_;
260 }
261
263
280 P1VTKFunction(const GV& gv, const V &v_, const std::string &s_,
281 int ncomps=1, int mycomp=0, VTK::Precision prec = VTK::Precision::float32)
282 : v( v_ ),
283 s( s_ ),
284 ncomps_(ncomps),
285 mycomp_(mycomp),
286 prec_(prec),
287 mapper( gv, mcmgVertexLayout() )
288 {
289 if (v.size()!=(unsigned int)(mapper.size()*ncomps_))
290 DUNE_THROW(IOError,"P1VTKFunction: size mismatch");
291 }
292
294 virtual ~P1VTKFunction() {}
295 };
296
298
299} // namespace Dune
300
301#endif // DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
@ mydimension
Dimensionality of the reference element of the entity.
Definition: entity.hh:113
Default exception class for I/O errors.
Definition: exceptions.hh:229
generic geometry implementation based on corner coordinates
Definition: multilineargeometry.hh:179
GlobalCoordinate global(const LocalCoordinate &local) const
evaluate the mapping
Definition: multilineargeometry.hh:282
Implementation class for a multiple codim and multiple geometry type mapper.
Definition: mcmgmapper.hh:200
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:268
size_type size() const
Return total number of entities in the entity set managed by the mapper.
Definition: mcmgmapper.hh:287
Index index(const EntityType &e) const
Map entity to starting index in array for dof block.
Definition: mcmgmapper.hh:254
Take a vector and interpret it as cell data for the VTKWriter.
Definition: function.hh:95
std::string name() const override
get name
Definition: function.hh:134
double evaluate(int, const Entity &e, const Dune::FieldVector< ctype, dim > &) const override
evaluate
Definition: function.hh:127
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:163
VTK::Precision precision() const override
get output precision for the field
Definition: function.hh:140
virtual ~P0VTKFunction()
destructor
Definition: function.hh:177
int ncomps() const override
return number of components
Definition: function.hh:121
Take a vector and interpret it as point data for the VTKWriter.
Definition: function.hh:203
std::string name() const override
get name
Definition: function.hh:251
virtual ~P1VTKFunction()
destructor
Definition: function.hh:294
VTK::Precision precision() const override
get output precision for the field
Definition: function.hh:257
double evaluate(int comp, const Entity &e, const Dune::FieldVector< ctype, dim > &xi) const override
evaluate
Definition: function.hh:235
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:280
int ncomps() const override
return number of components
Definition: function.hh:229
A base class for grid functions with any return type and dimension.
Definition: function.hh:40
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
virtual std::string name() const =0
get name
virtual VTK::Precision precision() const
get output precision for the field
Definition: function.hh:65
virtual int ncomps() const =0
virtual ~VTKFunction()
virtual destructor
Definition: function.hh:69
A few common exception classes.
Common stuff for the VTKWriter.
Precision
which precision to use when writing out data to vtk files
Definition: common.hh:319
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Grid::ctype ctype
type used for coordinates in grid
Definition: gridview.hh:124
@ dimension
The dimension of the grid.
Definition: gridview.hh:127
MCMGLayout mcmgElementLayout()
layout for elements (codim-0 entities)
Definition: mcmgmapper.hh:150
MCMGLayout mcmgVertexLayout()
layout for vertices (dim-0 entities)
Definition: mcmgmapper.hh:160
Mapper for multiple codim and multiple geometry types.
Dune namespace.
Definition: alignedallocator.hh:14
Static tag representing a codimension.
Definition: dimension.hh:22
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)