DUNE-FEM (unstable)

quadratureimp.hh
1#ifndef DUNE_FEM_QUADRATUREIMP_HH
2#define DUNE_FEM_QUADRATUREIMP_HH
3
4#include <cassert>
5#include <vector>
6
8
10
11#include <dune/fem/version.hh>
12#include <dune/fem/quadrature/idprovider.hh>
13
14namespace Dune
15{
16
17 namespace Fem
18 {
19
31 template< typename FieldImp, int dim >
33 {
34 public:
36 typedef FieldImp FieldType;
37
38 protected:
40
41 // this is used for interpolation quadratures to obtain element
42 // interpolation points for face quadratures
44 public:
47
49 enum { codimension = 0 };
50
52 static const int dimension = dim ;
53
54 protected:
55 // vector holding the coordinates for each point
56 mutable std::vector< CoordinateType > points_;
57
58 // identifier of the integration point list
59 const size_t id_;
60
76 explicit IntegrationPointListImp( size_t id )
77 : points_(), id_( id )
78 {}
79
80 public:
82
83 virtual ~IntegrationPointListImp () = default;
84
96 const CoordinateType &point ( size_t i ) const
97 {
98 assert( i < nop() );
99 return points_[ i ];
100 }
101
106 size_t nop () const
107 {
108 return points_.size();
109 }
110
122 size_t id () const
123 {
124 return id_;
125 }
126
137 virtual int order() const = 0;
138
148 virtual GeometryType geometryType () const = 0;
149
152 virtual std::vector< ElementCoordinateType > interpolationPoints(const int reqDim) const
153 {
154 return std::vector< ElementCoordinateType >();
155 }
156
160 virtual bool isFaceInterpolationQuadrature(const size_t numShapeFunctions ) const { return false; }
161
162 protected:
170 {
171 points_.push_back( point );
172 }
173
175 void setIntegrationPoints( std::vector< CoordinateType >&& points )
176 {
177 points_.swap( points );
178 }
179 };
180
181
182
193 template< typename FieldImp, int dim >
195 : public IntegrationPointListImp< FieldImp, dim >
196 {
197 public:
199 typedef FieldImp FieldType;
200
201 private:
204
205 public:
208
209 protected:
210 // vector holding weights of each integration point
211 mutable std::vector< FieldType > weights_;
212
228 explicit QuadratureImp( size_t id )
229 : BaseType( id ), weights_()
230 {}
231
232 public:
233 QuadratureImp ( const QuadratureImp& ) = delete;
234
235 virtual ~QuadratureImp () = default;
236
251 const FieldType &weight ( size_t i ) const
252 {
253 return weights_[ i ];
254 }
255
256 private:
257 // Disallow use of addIntegrationPoint for quadratures
258 void addIntegrationPoint ( const CoordinateType &point )
259 {
261 }
262
263 protected:
271 {
272 addIntegrationPoint( point );
273 weights_.push_back( weight );
274 }
275 };
276
277
278
279 // \brief Allows injection of arbitrary points as quadrature points.
280 // Useful to test some features of the quadrature framework in isolation
281 // and with known input data. Each TestQuadrature object gets its own
282 // unique id.
283 template <class ct, int dim>
284 class TestQuadrature : public QuadratureImp<ct, dim>
285 {
286 public:
287 typedef FieldVector<ct, dim> CoordinateType;
288
289 // dummy value
290 enum { maxOrder_ = 10 };
291
293 TestQuadrature(const GeometryType& geo, int order);
294
296 void newQuadraturePoint(const CoordinateType& c, ct weight);
297
299 virtual GeometryType geometryType() const { return geo_; }
300
302 virtual int order() const { return order_; }
303
305 static size_t maxOrder() { return maxOrder_; }
306
307 private:
308 GeometryType geo_;
309 int order_;
310 };
311
312 } // namespace Fem
313
314} // namespace Dune
315
316#include "quadratureimp_inline.hh"
317
318#endif // #ifndef DUNE_FEM_QUADRATUREIMP_HH
Generic implementation of an IntegrationPointList.
Definition: quadratureimp.hh:33
virtual std::vector< ElementCoordinateType > interpolationPoints(const int reqDim) const
returns list of element interpolation points for a given face quadrature
Definition: quadratureimp.hh:152
FieldImp FieldType
field type
Definition: quadratureimp.hh:36
void addIntegrationPoint(const CoordinateType &point)
Adds an integration point to the list.
Definition: quadratureimp.hh:169
virtual bool isFaceInterpolationQuadrature(const size_t numShapeFunctions) const
return true if quadrature is also a set of interpolation points for a given number of shape functions
Definition: quadratureimp.hh:160
size_t id() const
obtain the identifier of the integration point list
Definition: quadratureimp.hh:122
const CoordinateType & point(size_t i) const
obtain coordinates of i-th integration point
Definition: quadratureimp.hh:96
IntegrationPointListImp(size_t id)
Constructor.
Definition: quadratureimp.hh:76
static const int dimension
dimension of quadrature
Definition: quadratureimp.hh:52
void setIntegrationPoints(std::vector< CoordinateType > &&points)
Overwrites integration point list
Definition: quadratureimp.hh:175
size_t nop() const
obtain the number of integration points
Definition: quadratureimp.hh:106
FieldVector< FieldType, dim > CoordinateType
type of local coordinates
Definition: quadratureimp.hh:46
virtual GeometryType geometryType() const =0
obtain GeometryType for this integration point list
virtual int order() const =0
obtain order of the integration point list
Generic implementation of a Dune quadrature.
Definition: quadratureimp.hh:196
QuadratureImp(size_t id)
Constructor.
Definition: quadratureimp.hh:228
void addQuadraturePoint(const CoordinateType &point, const FieldType weight)
Adds a point-weight pair to the quadrature.
Definition: quadratureimp.hh:270
BaseType::CoordinateType CoordinateType
type of local coordinates
Definition: quadratureimp.hh:207
FieldImp FieldType
field type
Definition: quadratureimp.hh:199
const FieldType & weight(size_t i) const
obtain weight of i-th integration point
Definition: quadratureimp.hh:251
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
Implements a vector constructed from a given type representing a field and a compile-time given size.
Dune namespace.
Definition: alignedallocator.hh:13
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 (Jul 24, 22:29, 2024)