Dune Core Modules (2.6.0)

quadraturerules.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_GEOMETRY_QUADRATURERULES_HH
5#define DUNE_GEOMETRY_QUADRATURERULES_HH
6
7#include <algorithm>
8#include <iostream>
9#include <limits>
10#include <mutex>
11#include <utility>
12#include <vector>
13
17#include <dune/common/stdthread.hh>
19
20#include <dune/geometry/quadraturerules/nocopyvector.hh>
21#include <dune/geometry/type.hh>
23
29namespace Dune {
30
36
42 template<typename ct, int dim>
44 public:
46 enum { dimension = dim };
47
49 typedef ct Field;
50
53
55 QuadraturePoint (const Vector& x, ct w) : local(x)
56 {
57 weight_ = w;
58 }
59
61 const Vector& position () const
62 {
63 return local;
64 }
65
67 const ct &weight () const
68 {
69 return weight_;
70 }
71
72 protected:
74 ct weight_;
75 };
76
80 namespace QuadratureType {
81 enum Enum {
82 GaussLegendre = 0,
83
84 GaussJacobi_1_0 = 1,
85 GaussJacobi_2_0 = 2,
86
87 GaussLobatto = 4,
88 size
89 };
90 }
91
95 template<typename ct, int dim>
96 class QuadratureRule : public std::vector<QuadraturePoint<ct,dim> >
97 {
98 public:
104 QuadratureRule() : delivered_order(-1) {}
105
106 protected:
108 QuadratureRule(GeometryType t) : geometry_type(t), delivered_order(-1) {}
109
111 QuadratureRule(GeometryType t, int order) : geometry_type(t), delivered_order(order) {}
112 public:
114 enum { d=dim };
115
117 typedef ct CoordType;
118
120 virtual int order () const { return delivered_order; }
121
123 virtual GeometryType type () const { return geometry_type; }
124 virtual ~QuadratureRule(){}
125
128 typedef typename std::vector<QuadraturePoint<ct,dim> >::const_iterator iterator;
129
130 protected:
131 GeometryType geometry_type;
132 int delivered_order;
133 };
134
135 // Forward declaration of the factory class,
136 // needed internally by the QuadratureRules container class.
137 template<typename ctype, int dim> class QuadratureRuleFactory;
138
142 template<typename ctype, int dim>
144
149 static void initQuadratureRule(QuadratureRule *qr, QuadratureType::Enum qt,
150 const GeometryType &t, int p)
151 {
153 }
154
155 typedef NoCopyVector<std::pair<std::once_flag, QuadratureRule> >
156 QuadratureOrderVector; // indexed by quadrature order
159 static void initQuadratureOrderVector(QuadratureOrderVector *qov,
160 QuadratureType::Enum qt,
161 const GeometryType &t)
162 {
163 if(dim == 0)
164 // we only need one quadrature rule for points, not maxint
165 qov->resize(1);
166 else
168 }
169
170 typedef NoCopyVector<std::pair<std::once_flag, QuadratureOrderVector> >
171 GeometryTypeVector; // indexed by geometry type
174 static void initGeometryTypeVector(GeometryTypeVector *gtv)
175 {
176 gtv->resize(LocalGeometryTypeIndex::size(dim));
177 }
178
180 DUNE_EXPORT const QuadratureRule& _rule(const GeometryType& t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
181 {
182 assert(t.dim()==dim);
183
184 DUNE_ASSERT_CALL_ONCE();
185
186 static NoCopyVector<std::pair< // indexed by quadrature type
187 std::once_flag,
188 GeometryTypeVector
189 > > quadratureCache(QuadratureType::size);
190
191 auto & quadratureTypeLevel = quadratureCache[qt];
192 std::call_once(quadratureTypeLevel.first, initGeometryTypeVector,
193 &quadratureTypeLevel.second);
194
195 auto & geometryTypeLevel =
196 quadratureTypeLevel.second[LocalGeometryTypeIndex::index(t)];
197 std::call_once(geometryTypeLevel.first, initQuadratureOrderVector,
198 &geometryTypeLevel.second, qt, t);
199
200 // we only have one quadrature rule for points
201 auto & quadratureOrderLevel = geometryTypeLevel.second[dim == 0 ? 0 : p];
202 std::call_once(quadratureOrderLevel.first, initQuadratureRule,
203 &quadratureOrderLevel.second, qt, t, p);
204
205 return quadratureOrderLevel.second;
206 }
208 DUNE_EXPORT static QuadratureRules& instance()
209 {
210 static QuadratureRules instance;
211 return instance;
212 }
214 QuadratureRules () {}
215 public:
217 static unsigned
219 QuadratureType::Enum qt=QuadratureType::GaussLegendre)
220 {
222 }
223
225 static const QuadratureRule& rule(const GeometryType& t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
226 {
227 return instance()._rule(t,p,qt);
228 }
229
230#pragma GCC diagnostic push
231#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
233 static const QuadratureRule& rule(const GeometryType::BasicType t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
234 {
235 GeometryType gt(t,dim);
236 return instance()._rule(gt,p,qt);
237 }
238 };
239#pragma GCC diagnostic pop
240
241} // end namespace Dune
242
243#include "quadraturerules/pointquadrature.hh"
244
245namespace Dune {
246
248 template<typename ct, bool fundamental = std::numeric_limits<ct>::is_specialized>
250 template<typename ct>
251 struct GaussQuadratureInitHelper<ct, true> {
252 static void init(int p,
253 std::vector< FieldVector<ct, 1> > & _points,
254 std::vector< ct > & _weight,
255 int & delivered_order);
256 };
257 template<typename ct>
258 struct GaussQuadratureInitHelper<ct, false> {
259 static void init(int p,
260 std::vector< FieldVector<ct, 1> > & _points,
261 std::vector< ct > & _weight,
262 int & delivered_order);
263 };
264
266 template<typename ct>
268 public QuadratureRule<ct,1>
269 {
270 public:
271 // compile time parameters
272 enum { dim=1 };
273 enum { highest_order=61 };
274
276 private:
277 friend class QuadratureRuleFactory<ct,dim>;
280 {
282 std::vector< FieldVector<ct, dim> > _points;
283 std::vector< ct > _weight;
284
286 (p, _points, _weight, this->delivered_order);
287
288 assert(_points.size() == _weight.size());
289 for (size_t i = 0; i < _points.size(); i++)
290 this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
291 }
292 };
293
296
297} // namespace Dune
298
299#define DUNE_INCLUDING_IMPLEMENTATION
300#include "quadraturerules/gauss_imp.hh"
301
302namespace Dune {
303
305 template<typename ct,
306 bool fundamental = std::numeric_limits<ct>::is_specialized>
308 template<typename ct>
309 struct Jacobi1QuadratureInitHelper<ct, true> {
310 static void init(int p,
311 std::vector< FieldVector<ct, 1> > & _points,
312 std::vector< ct > & _weight,
313 int & delivered_order);
314 };
315 template<typename ct>
316 struct Jacobi1QuadratureInitHelper<ct, false> {
317 static void init(int p,
318 std::vector< FieldVector<ct, 1> > & _points,
319 std::vector< ct > & _weight,
320 int & delivered_order);
321 };
322
326 template<typename ct>
328 public QuadratureRule<ct,1>
329 {
330 public:
332 enum { dim=1 };
333
335 enum { highest_order=61 };
336
338 private:
339 friend class QuadratureRuleFactory<ct,dim>;
342 {
344 std::vector< FieldVector<ct, dim> > _points;
345 std::vector< ct > _weight;
346
347 int deliveredOrder_;
348
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++)
354 this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
355 }
356 };
357
358#ifndef DOXYGEN
361#endif // !DOXYGEN
362
363} // namespace Dune
364
365#define DUNE_INCLUDING_IMPLEMENTATION
366#include "quadraturerules/jacobi_1_0_imp.hh"
367
368namespace Dune {
369
371 template<typename ct,
372 bool fundamental = std::numeric_limits<ct>::is_specialized>
374 template<typename ct>
375 struct Jacobi2QuadratureInitHelper<ct, true> {
376 static void init(int p,
377 std::vector< FieldVector<ct, 1> > & _points,
378 std::vector< ct > & _weight,
379 int & delivered_order);
380 };
381 template<typename ct>
382 struct Jacobi2QuadratureInitHelper<ct, false> {
383 static void init(int p,
384 std::vector< FieldVector<ct, 1> > & _points,
385 std::vector< ct > & _weight,
386 int & delivered_order);
387 };
388
392 template<typename ct>
394 public QuadratureRule<ct,1>
395 {
396 public:
398 enum { dim=1 };
399
401 enum { highest_order=61 };
402
404 private:
405 friend class QuadratureRuleFactory<ct,dim>;
408 {
410 std::vector< FieldVector<ct, dim> > _points;
411 std::vector< ct > _weight;
412
413 int deliveredOrder_;
414
416 (p, _points, _weight, deliveredOrder_);
417
418 this->delivered_order = deliveredOrder_;
419 assert(_points.size() == _weight.size());
420 for (size_t i = 0; i < _points.size(); i++)
421 this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
422 }
423 };
424
425#ifndef DOXYGEN
428#endif // !DOXYGEN
429
430} // namespace Dune
431
432#define DUNE_INCLUDING_IMPLEMENTATION
433#include "quadraturerules/jacobi_2_0_imp.hh"
434
435namespace Dune {
436
438 template<typename ct,
439 bool fundamental = std::numeric_limits<ct>::is_specialized>
441 template<typename ct>
442 struct GaussLobattoQuadratureInitHelper<ct, true> {
443 static void init(int p,
444 std::vector< FieldVector<ct, 1> > & _points,
445 std::vector< ct > & _weight,
446 int & delivered_order);
447 };
448 template<typename ct>
449 struct GaussLobattoQuadratureInitHelper<ct, false> {
450 static void init(int p,
451 std::vector< FieldVector<ct, 1> > & _points,
452 std::vector< ct > & _weight,
453 int & delivered_order);
454 };
455
459 template<typename ct>
461 public QuadratureRule<ct,1>
462 {
463 public:
465 enum { dim=1 };
466
468 enum { highest_order=31 };
469
471 private:
472 friend class QuadratureRuleFactory<ct,dim>;
475 {
477 std::vector< FieldVector<ct, dim> > _points;
478 std::vector< ct > _weight;
479
480 int deliveredOrder_;
481
483 (p, _points, _weight, deliveredOrder_);
484
485 this->delivered_order = deliveredOrder_;
486 assert(_points.size() == _weight.size());
487 for (size_t i = 0; i < _points.size(); i++)
488 this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
489 }
490 };
491
492#ifndef DOXYGEN
495#endif // !DOXYGEN
496
497} // namespace Dune
498
499#define DUNE_INCLUDING_IMPLEMENTATION
500#include "quadraturerules/gausslobatto_imp.hh"
501
502#include "quadraturerules/tensorproductquadrature.hh"
503
504#include "quadraturerules/simplexquadrature.hh"
505
506namespace Dune {
507
508 /***********************************
509 * quadrature for Prism
510 **********************************/
511
513 template<int dim>
515
517 template<>
519 {
520 public:
521 enum { MAXP=6};
522 enum { highest_order=2 };
523
526 {
527 int m = 0;
528 O[m] = 0;
529
530 // polynom degree 0 ???
531 m = 6;
532 G[m][0][0] = 0.0;
533 G[m][0][1] = 0.0;
534 G[m][0][2] = 0.0;
535
536 G[m][1][0] = 1.0;
537 G[m][1][1] = 0.0;
538 G[m][1][2] = 0.0;
539
540 G[m][2][0] = 0.0;
541 G[m][2][1] = 1.0;
542 G[m][2][2] = 0.0;
543
544 G[m][3][0] = 0.0;
545 G[m][3][1] = 0.0;
546 G[m][3][2] = 1.0;
547
548 G[m][4][0] = 1.0;
549 G[m][4][1] = 0.0;
550 G[m][4][2] = 1.0;
551
552 G[m][5][0] = 0.0;
553 G[m][5][1] = 0.1;
554 G[m][5][2] = 1.0;
555
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;
562
563 O[m] = 0; // verify ????????
564
565
566 // polynom degree 2 ???
567 m = 6;
568 G[m][0][0] =0.66666666666666666 ;
569 G[m][0][1] =0.16666666666666666 ;
570 G[m][0][2] =0.211324865405187 ;
571
572 G[m][1][0] = 0.16666666666666666;
573 G[m][1][1] =0.66666666666666666 ;
574 G[m][1][2] = 0.211324865405187;
575
576 G[m][2][0] = 0.16666666666666666;
577 G[m][2][1] = 0.16666666666666666;
578 G[m][2][2] = 0.211324865405187;
579
580 G[m][3][0] = 0.66666666666666666;
581 G[m][3][1] = 0.16666666666666666;
582 G[m][3][2] = 0.788675134594813;
583
584 G[m][4][0] = 0.16666666666666666;
585 G[m][4][1] = 0.66666666666666666;
586 G[m][4][2] = 0.788675134594813;
587
588 G[m][5][0] = 0.16666666666666666;
589 G[m][5][1] = 0.16666666666666666;
590 G[m][5][2] = 0.788675134594813;
591
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;
598
599 O[m] = 2; // verify ????????
600
601 }
602
605 {
606 return G[m][i];
607 }
608
610 double weight (int m, int i)
611 {
612 return W[m][i];
613 }
614
616 int order (int m)
617 {
618 return O[m];
619 }
620
621 private:
622 FieldVector<double, 3> G[MAXP+1][MAXP]; //positions
623
624 double W[MAXP+1][MAXP]; // weights associated with points
625 int O[MAXP+1]; // order of the rule
626 };
627
628
632 template<int dim>
634 static PrismQuadraturePoints<3> prqp;
635 };
636
640 template<>
642 static PrismQuadraturePoints<3> prqp;
643 };
644
648 template<typename ct, int dim>
650
654 template<typename ct>
656 {
657 public:
658
660 enum { d = 3 };
661
663 enum { highest_order = 2 };
664
666 private:
667 friend class QuadratureRuleFactory<ct,d>;
669 {
670 int m=6;
671 this->delivered_order = PrismQuadraturePointsSingleton<3>::prqp.order(m);
672 for(int i=0; i<m; ++i)
673 {
674 FieldVector<ct,3> local;
675 for (int k=0; k<d; k++)
676 local[k] = PrismQuadraturePointsSingleton<3>::prqp.point(m,i)[k];
677 double weight =
679 // put in container
680 this->push_back(QuadraturePoint<ct,d>(local,weight));
681 }
682 }
683 };
684
691 template<typename ctype, int dim>
693 private:
694 friend class QuadratureRules<ctype, dim>;
695 static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
696 {
697 return TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
698 }
699 static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
700 {
701 return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
702 }
703 };
704
705 template<typename ctype>
707 private:
708 enum { dim = 0 };
709 friend class QuadratureRules<ctype, dim>;
710 static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
711 {
712 if (t.isVertex())
713 {
714 return std::numeric_limits<int>::max();
715 }
716 DUNE_THROW(Exception, "Unknown GeometryType");
717 }
718 static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
719 {
720 if (t.isVertex())
721 {
722 return PointQuadratureRule<ctype>();
723 }
724 DUNE_THROW(Exception, "Unknown GeometryType");
725 }
726 };
727
728 template<typename ctype>
729 class QuadratureRuleFactory<ctype, 1> {
730 private:
731 enum { dim = 1 };
732 friend class QuadratureRules<ctype, dim>;
733 static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
734 {
735 if (t.isLine())
736 {
737 switch (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;
746 default :
747 DUNE_THROW(Exception, "Unknown QuadratureType");
748 }
749 }
750 DUNE_THROW(Exception, "Unknown GeometryType");
751 }
752 static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
753 {
754 if (t.isLine())
755 {
756 switch (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);
765 default :
766 DUNE_THROW(Exception, "Unknown QuadratureType");
767 }
768 }
769 DUNE_THROW(Exception, "Unknown GeometryType");
770 }
771 };
772
773 template<typename ctype>
774 class QuadratureRuleFactory<ctype, 2> {
775 private:
776 enum { dim = 2 };
777 friend class QuadratureRules<ctype, dim>;
778 static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
779 {
780 unsigned order =
781 TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
782 if (t.isSimplex())
783 order = std::max
784 (order, unsigned(SimplexQuadratureRule<ctype,dim>::highest_order));
785 return order;
786 }
787 static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
788 {
789 if (t.isSimplex()
790 && qt == QuadratureType::GaussLegendre
791 && p <= SimplexQuadratureRule<ctype,dim>::highest_order)
792 {
793 return SimplexQuadratureRule<ctype,dim>(p);
794 }
795 return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
796 }
797 };
798
799 template<typename ctype>
800 class QuadratureRuleFactory<ctype, 3> {
801 private:
802 enum { dim = 3 };
803 friend class QuadratureRules<ctype, dim>;
804 static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
805 {
806 unsigned order =
807 TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
808 if (t.isSimplex())
809 order = std::max
810 (order, unsigned(SimplexQuadratureRule<ctype,dim>::highest_order));
811 if (t.isPrism())
812 order = std::max
813 (order, unsigned(PrismQuadratureRule<ctype,dim>::highest_order));
814 return order;
815 }
816 static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
817 {
818 if (t.isSimplex()
819 && qt == QuadratureType::GaussLegendre
820 && p <= SimplexQuadratureRule<ctype,dim>::highest_order)
821 {
822 return SimplexQuadratureRule<ctype,dim>(p);
823 }
824 if (t.isPrism()
825 && qt == QuadratureType::GaussLegendre
826 && p <= PrismQuadratureRule<ctype,dim>::highest_order)
827 {
828 return PrismQuadratureRule<ctype,dim>(p);
829 }
830 return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
831 }
832 };
833
834} // end namespace
835
836#endif // DUNE_GEOMETRY_QUADRATURERULES_HH
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 &gt)
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)