DUNE PDELab (git)

function.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5
6#ifndef DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
7#define DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
8
9#include <string>
10
13
14#include <dune/geometry/type.hh>
15#include <dune/geometry/referenceelements.hh>
16#include <dune/geometry/multilineargeometry.hh>
17
20
26namespace Dune
27{
30
32 //
33 // Base VTKFunction
34 //
35
40 template< class GridView >
42 {
43 public:
44 typedef typename GridView::ctype ctype;
45 constexpr static int dim = GridView::dimension;
46 typedef typename GridView::template Codim< 0 >::Entity Entity;
47
50 virtual int ncomps () const = 0;
51
53
60 virtual double evaluate (int comp, const Entity& e,
61 const Dune::FieldVector<ctype,dim>& xi) const = 0;
62
64 virtual std::string name () const = 0;
65
67 virtual VTK::Precision precision() const
68 { return VTK::Precision::float32; }
69
71 virtual ~VTKFunction () {}
72 };
73
75 //
76 // P0VTKFunction
77 //
78
80
94 template<typename GV, typename V>
96 : public VTKFunction< GV >
97 {
99 typedef VTKFunction< GV > Base;
102
104 const V& v;
106 std::string s;
108 int ncomps_;
111 int mycomp_;
113 VTK::Precision prec_;
115 Mapper mapper;
116
117 public:
118 typedef typename Base::Entity Entity;
119 typedef typename Base::ctype ctype;
120 using Base::dim;
121
123 int ncomps () const override
124 {
125 return 1;
126 }
127
129 double evaluate (int, const Entity& e,
130 const Dune::FieldVector<ctype,dim>&) const override
131 {
132 return v[mapper.index(e)*ncomps_+mycomp_];
133 }
134
136 std::string name () const override
137 {
138 return s;
139 }
140
142 VTK::Precision precision() const override
143 {
144 return prec_;
145 }
146
148
165 P0VTKFunction(const GV &gv, const V &v_, const std::string &s_,
166 int ncomps=1, int mycomp=0, VTK::Precision prec = VTK::Precision::float32)
167 : v( v_ ),
168 s( s_ ),
169 ncomps_(ncomps),
170 mycomp_(mycomp),
171 prec_(prec),
172 mapper( gv, mcmgElementLayout() )
173 {
174 if (v.size()!=(unsigned int)(mapper.size()*ncomps_))
175 DUNE_THROW(IOError, "P0VTKFunction: size mismatch");
176 }
177
179 virtual ~P0VTKFunction() {}
180 };
181
183 //
184 // P1VTKFunction
185 //
186
188
202 template<typename GV, typename V>
204 : public VTKFunction< GV >
205 {
207 typedef VTKFunction< GV > Base;
210
212 const V& v;
214 std::string s;
216 int ncomps_;
219 int mycomp_;
221 VTK::Precision prec_;
223 Mapper mapper;
224
225 public:
226 typedef typename Base::Entity Entity;
227 typedef typename Base::ctype ctype;
228 using Base::dim;
229
231 int ncomps () const override
232 {
233 return 1;
234 }
235
237 double evaluate ([[maybe_unused]] int comp, const Entity& e,
238 const Dune::FieldVector<ctype,dim>& xi) const override
239 {
240 const unsigned int myDim = Entity::mydimension;
241 const unsigned int nVertices = e.subEntities(dim);
242
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_];
246
247 // (Ab)use the MultiLinearGeometry class to do multi-linear interpolation between scalars
248 const MultiLinearGeometry<ctype,dim,1> interpolation(e.type(), cornerValues);
249 return interpolation.global(xi);
250 }
251
253 std::string name () const override
254 {
255 return s;
256 }
257
259 VTK::Precision precision() const override
260 {
261 return prec_;
262 }
263
265
282 P1VTKFunction(const GV& gv, const V &v_, const std::string &s_,
283 int ncomps=1, int mycomp=0, VTK::Precision prec = VTK::Precision::float32)
284 : v( v_ ),
285 s( s_ ),
286 ncomps_(ncomps),
287 mycomp_(mycomp),
288 prec_(prec),
289 mapper( gv, mcmgVertexLayout() )
290 {
291 if (v.size()!=(unsigned int)(mapper.size()*ncomps_))
292 DUNE_THROW(IOError,"P1VTKFunction: size mismatch");
293 }
294
296 virtual ~P1VTKFunction() {}
297 };
298
300
301} // namespace Dune
302
303#endif // DUNE_GRID_IO_FILE_VTK_FUNCTION_HH
static constexpr int mydimension
Dimensionality of the reference element of the entity.
Definition: entity.hh:112
Default exception class for I/O errors.
Definition: exceptions.hh:231
generic geometry implementation based on corner coordinates
Definition: multilineargeometry.hh:181
GlobalCoordinate global(const LocalCoordinate &local) const
evaluate the mapping
Definition: multilineargeometry.hh:290
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
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
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
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
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
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
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
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:271
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:218
static constexpr int dimension
The dimension of the grid.
Definition: gridview.hh:134
Grid::ctype ctype
type used for coordinates in grid
Definition: gridview.hh:131
MCMGLayout mcmgElementLayout()
layout for elements (codim-0 entities)
Definition: mcmgmapper.hh:97
MCMGLayout mcmgVertexLayout()
layout for vertices (dim-0 entities)
Definition: mcmgmapper.hh:107
Mapper for multiple codim and multiple geometry types.
Dune namespace.
Definition: alignedallocator.hh:13
Static tag representing a codimension.
Definition: dimension.hh:24
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 24, 23:30, 2024)