4#ifndef DUNE_GEOMETRY_QUADRATURERULES_HH
5#define DUNE_GEOMETRY_QUADRATURERULES_HH
17#include <dune/common/stdthread.hh>
20#include <dune/geometry/quadraturerules/nocopyvector.hh>
42 template<
typename ct,
int dim>
46 enum { dimension = dim };
80 namespace QuadratureType {
95 template<
typename ct,
int dim>
120 virtual int order ()
const {
return delivered_order; }
128 typedef typename std::vector<QuadraturePoint<ct,dim> >::const_iterator
iterator;
142 template<
typename ctype,
int dim>
149 static void initQuadratureRule(
QuadratureRule *qr, QuadratureType::Enum qt,
155 typedef NoCopyVector<std::pair<std::once_flag, QuadratureRule> >
156 QuadratureOrderVector;
159 static void initQuadratureOrderVector(QuadratureOrderVector *qov,
160 QuadratureType::Enum qt,
170 typedef NoCopyVector<std::pair<std::once_flag, QuadratureOrderVector> >
174 static void initGeometryTypeVector(GeometryTypeVector *gtv)
182 assert(t.
dim()==dim);
184 DUNE_ASSERT_CALL_ONCE();
186 static NoCopyVector<std::pair<
189 > > quadratureCache(QuadratureType::size);
191 auto & quadratureTypeLevel = quadratureCache[qt];
192 std::call_once(quadratureTypeLevel.first, initGeometryTypeVector,
193 &quadratureTypeLevel.second);
195 auto & geometryTypeLevel =
197 std::call_once(geometryTypeLevel.first, initQuadratureOrderVector,
198 &geometryTypeLevel.second, qt, t);
201 auto & quadratureOrderLevel = geometryTypeLevel.second[dim == 0 ? 0 : p];
202 std::call_once(quadratureOrderLevel.first, initQuadratureRule,
203 &quadratureOrderLevel.second, qt, t, p);
205 return quadratureOrderLevel.second;
219 QuadratureType::Enum qt=QuadratureType::GaussLegendre)
227 return instance()._rule(t,p,qt);
230#pragma GCC diagnostic push
231#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
236 return instance()._rule(
gt,p,qt);
239#pragma GCC diagnostic pop
243#include "quadraturerules/pointquadrature.hh"
248 template<typename ct, bool fundamental = std::numeric_limits<ct>::is_specialized>
250 template<
typename ct>
252 static void init(
int p,
254 std::vector< ct > & _weight,
255 int & delivered_order);
257 template<
typename ct>
259 static void init(
int p,
261 std::vector< ct > & _weight,
262 int & delivered_order);
266 template<
typename ct>
273 enum { highest_order=61 };
282 std::vector< FieldVector<ct, dim> > _points;
283 std::vector< ct > _weight;
286 (p, _points, _weight, this->delivered_order);
288 assert(_points.size() == _weight.size());
289 for (
size_t i = 0; i < _points.size(); i++)
299#define DUNE_INCLUDING_IMPLEMENTATION
300#include "quadraturerules/gauss_imp.hh"
305 template<
typename ct,
306 bool fundamental = std::numeric_limits<ct>::is_specialized>
308 template<
typename ct>
310 static void init(
int p,
312 std::vector< ct > & _weight,
313 int & delivered_order);
315 template<
typename ct>
317 static void init(
int p,
319 std::vector< ct > & _weight,
320 int & delivered_order);
326 template<
typename ct>
335 enum { highest_order=61 };
344 std::vector< FieldVector<ct, dim> > _points;
345 std::vector< ct > _weight;
350 (p, _points, _weight, deliveredOrder_);
351 this->delivered_order = deliveredOrder_;
352 assert(_points.size() == _weight.size());
353 for (
size_t i = 0; i < _points.size(); i++)
365#define DUNE_INCLUDING_IMPLEMENTATION
366#include "quadraturerules/jacobi_1_0_imp.hh"
371 template<
typename ct,
372 bool fundamental = std::numeric_limits<ct>::is_specialized>
374 template<
typename ct>
376 static void init(
int p,
378 std::vector< ct > & _weight,
379 int & delivered_order);
381 template<
typename ct>
383 static void init(
int p,
385 std::vector< ct > & _weight,
386 int & delivered_order);
392 template<
typename ct>
401 enum { highest_order=61 };
410 std::vector< FieldVector<ct, dim> > _points;
411 std::vector< ct > _weight;
416 (p, _points, _weight, deliveredOrder_);
418 this->delivered_order = deliveredOrder_;
419 assert(_points.size() == _weight.size());
420 for (
size_t i = 0; i < _points.size(); i++)
432#define DUNE_INCLUDING_IMPLEMENTATION
433#include "quadraturerules/jacobi_2_0_imp.hh"
438 template<
typename ct,
439 bool fundamental = std::numeric_limits<ct>::is_specialized>
441 template<
typename ct>
443 static void init(
int p,
445 std::vector< ct > & _weight,
446 int & delivered_order);
448 template<
typename ct>
450 static void init(
int p,
452 std::vector< ct > & _weight,
453 int & delivered_order);
459 template<
typename ct>
468 enum { highest_order=31 };
477 std::vector< FieldVector<ct, dim> > _points;
478 std::vector< ct > _weight;
483 (p, _points, _weight, deliveredOrder_);
485 this->delivered_order = deliveredOrder_;
486 assert(_points.size() == _weight.size());
487 for (
size_t i = 0; i < _points.size(); i++)
499#define DUNE_INCLUDING_IMPLEMENTATION
500#include "quadraturerules/gausslobatto_imp.hh"
502#include "quadraturerules/tensorproductquadrature.hh"
504#include "quadraturerules/simplexquadrature.hh"
522 enum { highest_order=2 };
556 W[m][0] = 0.16666666666666666 / 2.0;
557 W[m][1] = 0.16666666666666666 / 2.0;
558 W[m][2] = 0.16666666666666666 / 2.0;
559 W[m][3] = 0.16666666666666666 / 2.0;
560 W[m][4] = 0.16666666666666666 / 2.0;
561 W[m][5] = 0.16666666666666666 / 2.0;
568 G[m][0][0] =0.66666666666666666 ;
569 G[m][0][1] =0.16666666666666666 ;
570 G[m][0][2] =0.211324865405187 ;
572 G[m][1][0] = 0.16666666666666666;
573 G[m][1][1] =0.66666666666666666 ;
574 G[m][1][2] = 0.211324865405187;
576 G[m][2][0] = 0.16666666666666666;
577 G[m][2][1] = 0.16666666666666666;
578 G[m][2][2] = 0.211324865405187;
580 G[m][3][0] = 0.66666666666666666;
581 G[m][3][1] = 0.16666666666666666;
582 G[m][3][2] = 0.788675134594813;
584 G[m][4][0] = 0.16666666666666666;
585 G[m][4][1] = 0.66666666666666666;
586 G[m][4][2] = 0.788675134594813;
588 G[m][5][0] = 0.16666666666666666;
589 G[m][5][1] = 0.16666666666666666;
590 G[m][5][2] = 0.788675134594813;
592 W[m][0] = 0.16666666666666666 / 2.0;
593 W[m][1] = 0.16666666666666666 / 2.0;
594 W[m][2] = 0.16666666666666666 / 2.0;
595 W[m][3] = 0.16666666666666666 / 2.0;
596 W[m][4] = 0.16666666666666666 / 2.0;
597 W[m][5] = 0.16666666666666666 / 2.0;
624 double W[MAXP+1][MAXP];
648 template<
typename ct,
int dim>
654 template<
typename ct>
663 enum { highest_order = 2 };
672 for(
int i=0; i<m; ++i)
675 for (
int k=0; k<d; k++)
691 template<
typename ctype,
int dim>
695 static unsigned maxOrder(
const GeometryType &t, QuadratureType::Enum qt)
697 return TensorProductQuadratureRule<ctype,dim>::maxOrder(t.
id(), qt);
701 return TensorProductQuadratureRule<ctype,dim>(t.
id(), p, qt);
705 template<
typename ctype>
709 friend class QuadratureRules<ctype, dim>;
710 static unsigned maxOrder(
const GeometryType &t, QuadratureType::Enum qt)
714 return std::numeric_limits<int>::max();
716 DUNE_THROW(Exception,
"Unknown GeometryType");
718 static QuadratureRule<ctype, dim> rule(
const GeometryType& t,
int p, QuadratureType::Enum qt)
722 return PointQuadratureRule<ctype>();
724 DUNE_THROW(Exception,
"Unknown GeometryType");
728 template<
typename ctype>
729 class QuadratureRuleFactory<ctype, 1> {
732 friend class QuadratureRules<ctype, dim>;
733 static unsigned maxOrder(
const GeometryType &t, QuadratureType::Enum qt)
738 case QuadratureType::GaussLegendre :
739 return GaussQuadratureRule1D<ctype>::highest_order;
740 case QuadratureType::GaussJacobi_1_0 :
741 return Jacobi1QuadratureRule1D<ctype>::highest_order;
742 case QuadratureType::GaussJacobi_2_0 :
743 return Jacobi2QuadratureRule1D<ctype>::highest_order;
744 case QuadratureType::GaussLobatto :
745 return GaussLobattoQuadratureRule1D<ctype>::highest_order;
747 DUNE_THROW(Exception,
"Unknown QuadratureType");
750 DUNE_THROW(Exception,
"Unknown GeometryType");
752 static QuadratureRule<ctype, dim> rule(
const GeometryType& t,
int p, QuadratureType::Enum qt)
757 case QuadratureType::GaussLegendre :
758 return GaussQuadratureRule1D<ctype>(p);
759 case QuadratureType::GaussJacobi_1_0 :
760 return Jacobi1QuadratureRule1D<ctype>(p);
761 case QuadratureType::GaussJacobi_2_0 :
762 return Jacobi2QuadratureRule1D<ctype>(p);
763 case QuadratureType::GaussLobatto :
764 return GaussLobattoQuadratureRule1D<ctype>(p);
766 DUNE_THROW(Exception,
"Unknown QuadratureType");
769 DUNE_THROW(Exception,
"Unknown GeometryType");
773 template<
typename ctype>
774 class QuadratureRuleFactory<ctype, 2> {
777 friend class QuadratureRules<ctype, dim>;
778 static unsigned maxOrder(
const GeometryType &t, QuadratureType::Enum qt)
781 TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
784 (order,
unsigned(SimplexQuadratureRule<ctype,dim>::highest_order));
787 static QuadratureRule<ctype, dim> rule(
const GeometryType& t,
int p, QuadratureType::Enum qt)
790 && qt == QuadratureType::GaussLegendre
791 && p <= SimplexQuadratureRule<ctype,dim>::highest_order)
793 return SimplexQuadratureRule<ctype,dim>(p);
795 return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
799 template<
typename ctype>
800 class QuadratureRuleFactory<ctype, 3> {
803 friend class QuadratureRules<ctype, dim>;
804 static unsigned maxOrder(
const GeometryType &t, QuadratureType::Enum qt)
807 TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
810 (order,
unsigned(SimplexQuadratureRule<ctype,dim>::highest_order));
813 (order,
unsigned(PrismQuadratureRule<ctype,dim>::highest_order));
816 static QuadratureRule<ctype, dim> rule(
const GeometryType& t,
int p, QuadratureType::Enum qt)
819 && qt == QuadratureType::GaussLegendre
820 && p <= SimplexQuadratureRule<ctype,dim>::highest_order)
822 return SimplexQuadratureRule<ctype,dim>(p);
825 && qt == QuadratureType::GaussLegendre
826 && p <= PrismQuadratureRule<ctype,dim>::highest_order)
828 return PrismQuadratureRule<ctype,dim>(p);
830 return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
Jacobi-Gauss quadrature for alpha=2, beta=0.
Definition: quadraturerules.hh:462
Gauss quadrature rule in 1D.
Definition: quadraturerules.hh:269
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
constexpr unsigned int dim() const
Return dimension of the type.
Definition: type.hh:572
BasicType
Each entity can be tagged by one of these basic types plus its space dimension.
Definition: type.hh:283
constexpr unsigned int id() const
Return the topology id of the type.
Definition: type.hh:577
Jacobi-Gauss quadrature for alpha=1, beta=0.
Definition: quadraturerules.hh:329
Jacobi-Gauss quadrature for alpha=2, beta=0.
Definition: quadraturerules.hh:395
static constexpr std::size_t size(std::size_t dim)
Compute total number of geometry types for the given dimension.
Definition: typeindex.hh:56
static constexpr std::size_t index(const GeometryType >)
Compute the index for the given geometry type within its dimension.
Definition: typeindex.hh:68
Default exception for dummy implementations.
Definition: exceptions.hh:261
Definition: quadraturerules.hh:519
double weight(int m, int i)
Definition: quadraturerules.hh:610
int order(int m)
Definition: quadraturerules.hh:616
PrismQuadraturePoints()
initialize quadrature points on the interval for all orders
Definition: quadraturerules.hh:525
FieldVector< double, 3 > point(int m, int i)
Definition: quadraturerules.hh:604
Definition: quadraturerules.hh:514
Quadrature rules for prisms.
Definition: quadraturerules.hh:656
Quadrature rules for prisms.
Definition: quadraturerules.hh:649
Exception thrown if a desired QuadratureRule is not available, because the requested order is to high...
Definition: quadraturerules.hh:35
Single evaluation point in a quadrature rule.
Definition: quadraturerules.hh:43
const Vector & position() const
return local coordinates of integration point i
Definition: quadraturerules.hh:61
Dune::FieldVector< ct, dim > Vector
Type used for the position of a quadrature point.
Definition: quadraturerules.hh:52
ct Field
Number type used for coordinates and quadrature weights.
Definition: quadraturerules.hh:49
const ct & weight() const
return weight associated with integration point i
Definition: quadraturerules.hh:67
QuadraturePoint(const Vector &x, ct w)
set up quadrature of given order in d dimensions
Definition: quadraturerules.hh:55
Factory class for creation of quadrature rules, depending on GeometryType, order and QuadratureType.
Definition: quadraturerules.hh:692
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:97
virtual GeometryType type() const
return type of element
Definition: quadraturerules.hh:123
QuadratureRule(GeometryType t, int order)
Constructor for a given geometry type and a given quadrature order.
Definition: quadraturerules.hh:111
ct CoordType
The type used for coordinates.
Definition: quadraturerules.hh:117
QuadratureRule()
Default constructor.
Definition: quadraturerules.hh:104
virtual int order() const
return order
Definition: quadraturerules.hh:120
QuadratureRule(GeometryType t)
Constructor for a given geometry type. Leaves the quadrature order invalid
Definition: quadraturerules.hh:108
std::vector< QuadraturePoint< ct, dim > >::const_iterator iterator
Definition: quadraturerules.hh:128
A container for all quadrature rules of dimension dim
Definition: quadraturerules.hh:143
static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
maximum quadrature order for given geometry type and quadrature type
Definition: quadraturerules.hh:218
static const QuadratureRule & rule(const GeometryType::BasicType t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:233
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:225
A few common exception classes.
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
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:147
constexpr GeometryType line
GeometryType representing a line.
Definition: type.hh:733
constexpr GeometryType prism
GeometryType representing a 3D prism.
Definition: type.hh:763
Dune namespace.
Definition: alignedallocator.hh:10
Standard Dune debug streams.
Definition: quadraturerules.hh:440
Definition: quadraturerules.hh:249
Definition: quadraturerules.hh:307
Definition: quadraturerules.hh:373
Singleton holding the Prism Quadrature points.
Definition: quadraturerules.hh:641
Singleton holding the Prism Quadrature points.
Definition: quadraturerules.hh:633
A unique label for each type of element that can occur in a grid.
Helper classes to provide indices for geometrytypes for use in a vector.
Definition of macros controlling symbol visibility at the ABI level.
#define DUNE_EXPORT
Export a symbol as part of the public ABI.
Definition: visibility.hh:18