Dune Core Modules (2.4.1)

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 
14 #include <dune/common/fvector.hh>
17 #include <dune/common/stdthread.hh>
19 
20 #include <dune/geometry/quadraturerules/nocopyvector.hh>
21 #include <dune/geometry/type.hh>
23 
29 namespace 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
167  qov->resize(QuadratureRuleFactory<ctype,dim>::maxOrder(t,qt)+1);
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  }
230  static const QuadratureRule& rule(const GeometryType::BasicType t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
231  {
232  GeometryType gt(t,dim);
233  return instance()._rule(gt,p,qt);
234  }
235  };
236 
237 } // end namespace Dune
238 
239 #include "quadraturerules/pointquadrature.hh"
240 
241 namespace Dune {
242 
244  template<typename ct, bool fundamental = std::numeric_limits<ct>::is_specialized>
246  template<typename ct>
247  struct GaussQuadratureInitHelper<ct, true> {
248  static void init(int p,
249  std::vector< FieldVector<ct, 1> > & _points,
250  std::vector< ct > & _weight,
251  int & delivered_order);
252  };
253  template<typename ct>
254  struct GaussQuadratureInitHelper<ct, false> {
255  static void init(int p,
256  std::vector< FieldVector<ct, 1> > & _points,
257  std::vector< ct > & _weight,
258  int & delivered_order);
259  };
260 
262  template<typename ct>
264  public QuadratureRule<ct,1>
265  {
266  public:
267  // compile time parameters
268  enum { dim=1 };
269  enum { highest_order=61 };
270 
272  private:
273  friend class QuadratureRuleFactory<ct,dim>;
274  GaussQuadratureRule1D (int p)
276  {
278  std::vector< FieldVector<ct, dim> > _points;
279  std::vector< ct > _weight;
280 
282  (p, _points, _weight, this->delivered_order);
283 
284  assert(_points.size() == _weight.size());
285  for (size_t i = 0; i < _points.size(); i++)
286  this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
287  }
288  };
289 
292 
293 } // namespace Dune
294 
295 #define DUNE_INCLUDING_IMPLEMENTATION
296 #include "quadraturerules/gauss_imp.hh"
297 
298 namespace Dune {
299 
301  template<typename ct,
302  bool fundamental = std::numeric_limits<ct>::is_specialized>
304  template<typename ct>
305  struct Jacobi1QuadratureInitHelper<ct, true> {
306  static void init(int p,
307  std::vector< FieldVector<ct, 1> > & _points,
308  std::vector< ct > & _weight,
309  int & delivered_order);
310  };
311  template<typename ct>
312  struct Jacobi1QuadratureInitHelper<ct, false> {
313  static void init(int p,
314  std::vector< FieldVector<ct, 1> > & _points,
315  std::vector< ct > & _weight,
316  int & delivered_order);
317  };
318 
322  template<typename ct>
324  public QuadratureRule<ct,1>
325  {
326  public:
328  enum { dim=1 };
329 
331  enum { highest_order=61 };
332 
334  private:
335  friend class QuadratureRuleFactory<ct,dim>;
337  : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
338  {
340  std::vector< FieldVector<ct, dim> > _points;
341  std::vector< ct > _weight;
342 
343  int delivered_order;
344 
346  (p, _points, _weight, delivered_order);
347  this->delivered_order = delivered_order;
348  assert(_points.size() == _weight.size());
349  for (size_t i = 0; i < _points.size(); i++)
350  this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
351  }
352  };
353 
354 #ifndef DOXYGEN
357 #endif // !DOXYGEN
358 
359 } // namespace Dune
360 
361 #define DUNE_INCLUDING_IMPLEMENTATION
362 #include "quadraturerules/jacobi_1_0_imp.hh"
363 
364 namespace Dune {
365 
367  template<typename ct,
368  bool fundamental = std::numeric_limits<ct>::is_specialized>
370  template<typename ct>
371  struct Jacobi2QuadratureInitHelper<ct, true> {
372  static void init(int p,
373  std::vector< FieldVector<ct, 1> > & _points,
374  std::vector< ct > & _weight,
375  int & delivered_order);
376  };
377  template<typename ct>
378  struct Jacobi2QuadratureInitHelper<ct, false> {
379  static void init(int p,
380  std::vector< FieldVector<ct, 1> > & _points,
381  std::vector< ct > & _weight,
382  int & delivered_order);
383  };
384 
388  template<typename ct>
390  public QuadratureRule<ct,1>
391  {
392  public:
394  enum { dim=1 };
395 
397  enum { highest_order=61 };
398 
400  private:
401  friend class QuadratureRuleFactory<ct,dim>;
403  : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
404  {
406  std::vector< FieldVector<ct, dim> > _points;
407  std::vector< ct > _weight;
408 
409  int delivered_order;
410 
412  (p, _points, _weight, delivered_order);
413 
414  this->delivered_order = delivered_order;
415  assert(_points.size() == _weight.size());
416  for (size_t i = 0; i < _points.size(); i++)
417  this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
418  }
419  };
420 
421 #ifndef DOXYGEN
424 #endif // !DOXYGEN
425 
426 } // namespace Dune
427 
428 #define DUNE_INCLUDING_IMPLEMENTATION
429 #include "quadraturerules/jacobi_2_0_imp.hh"
430 
431 namespace Dune {
432 
434  template<typename ct,
435  bool fundamental = std::numeric_limits<ct>::is_specialized>
437  template<typename ct>
438  struct GaussLobattoQuadratureInitHelper<ct, true> {
439  static void init(int p,
440  std::vector< FieldVector<ct, 1> > & _points,
441  std::vector< ct > & _weight,
442  int & delivered_order);
443  };
444  template<typename ct>
445  struct GaussLobattoQuadratureInitHelper<ct, false> {
446  static void init(int p,
447  std::vector< FieldVector<ct, 1> > & _points,
448  std::vector< ct > & _weight,
449  int & delivered_order);
450  };
451 
455  template<typename ct>
457  public QuadratureRule<ct,1>
458  {
459  public:
461  enum { dim=1 };
462 
464  enum { highest_order=31 };
465 
467  private:
468  friend class QuadratureRuleFactory<ct,dim>;
470  : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
471  {
473  std::vector< FieldVector<ct, dim> > _points;
474  std::vector< ct > _weight;
475 
476  int delivered_order;
477 
479  (p, _points, _weight, delivered_order);
480 
481  this->delivered_order = delivered_order;
482  assert(_points.size() == _weight.size());
483  for (size_t i = 0; i < _points.size(); i++)
484  this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
485  }
486  };
487 
488 #ifndef DOXYGEN
491 #endif // !DOXYGEN
492 
493 } // namespace Dune
494 
495 #define DUNE_INCLUDING_IMPLEMENTATION
496 #include "quadraturerules/gausslobatto_imp.hh"
497 
498 #include "quadraturerules/tensorproductquadrature.hh"
499 
500 #include "quadraturerules/simplexquadrature.hh"
501 
502 namespace Dune {
503 
504  /***********************************
505  * quadrature for Prism
506  **********************************/
507 
509  template<int dim>
511 
513  template<>
515  {
516  public:
517  enum { MAXP=6};
518  enum { highest_order=2 };
519 
522  {
523  int m = 0;
524  O[m] = 0;
525 
526  // polynom degree 0 ???
527  m = 6;
528  G[m][0][0] = 0.0;
529  G[m][0][1] = 0.0;
530  G[m][0][2] = 0.0;
531 
532  G[m][1][0] = 1.0;
533  G[m][1][1] = 0.0;
534  G[m][1][2] = 0.0;
535 
536  G[m][2][0] = 0.0;
537  G[m][2][1] = 1.0;
538  G[m][2][2] = 0.0;
539 
540  G[m][3][0] = 0.0;
541  G[m][3][1] = 0.0;
542  G[m][3][2] = 1.0;
543 
544  G[m][4][0] = 1.0;
545  G[m][4][1] = 0.0;
546  G[m][4][2] = 1.0;
547 
548  G[m][5][0] = 0.0;
549  G[m][5][1] = 0.1;
550  G[m][5][2] = 1.0;
551 
552  W[m][0] = 0.16666666666666666 / 2.0;
553  W[m][1] = 0.16666666666666666 / 2.0;
554  W[m][2] = 0.16666666666666666 / 2.0;
555  W[m][3] = 0.16666666666666666 / 2.0;
556  W[m][4] = 0.16666666666666666 / 2.0;
557  W[m][5] = 0.16666666666666666 / 2.0;
558 
559  O[m] = 0; // verify ????????
560 
561 
562  // polynom degree 2 ???
563  m = 6;
564  G[m][0][0] =0.66666666666666666 ;
565  G[m][0][1] =0.16666666666666666 ;
566  G[m][0][2] =0.211324865405187 ;
567 
568  G[m][1][0] = 0.16666666666666666;
569  G[m][1][1] =0.66666666666666666 ;
570  G[m][1][2] = 0.211324865405187;
571 
572  G[m][2][0] = 0.16666666666666666;
573  G[m][2][1] = 0.16666666666666666;
574  G[m][2][2] = 0.211324865405187;
575 
576  G[m][3][0] = 0.66666666666666666;
577  G[m][3][1] = 0.16666666666666666;
578  G[m][3][2] = 0.788675134594813;
579 
580  G[m][4][0] = 0.16666666666666666;
581  G[m][4][1] = 0.66666666666666666;
582  G[m][4][2] = 0.788675134594813;
583 
584  G[m][5][0] = 0.16666666666666666;
585  G[m][5][1] = 0.16666666666666666;
586  G[m][5][2] = 0.788675134594813;
587 
588  W[m][0] = 0.16666666666666666 / 2.0;
589  W[m][1] = 0.16666666666666666 / 2.0;
590  W[m][2] = 0.16666666666666666 / 2.0;
591  W[m][3] = 0.16666666666666666 / 2.0;
592  W[m][4] = 0.16666666666666666 / 2.0;
593  W[m][5] = 0.16666666666666666 / 2.0;
594 
595  O[m] = 2; // verify ????????
596 
597  }
598 
601  {
602  return G[m][i];
603  }
604 
606  double weight (int m, int i)
607  {
608  return W[m][i];
609  }
610 
612  int order (int m)
613  {
614  return O[m];
615  }
616 
617  private:
618  FieldVector<double, 3> G[MAXP+1][MAXP]; //positions
619 
620  double W[MAXP+1][MAXP]; // weights associated with points
621  int O[MAXP+1]; // order of the rule
622  };
623 
624 
628  template<int dim>
630  static PrismQuadraturePoints<3> prqp;
631  };
632 
636  template<>
638  static PrismQuadraturePoints<3> prqp;
639  };
640 
644  template<typename ct, int dim>
646 
650  template<typename ct>
652  {
653  public:
654 
656  enum { d = 3 };
657 
659  enum { highest_order = 2 };
660 
662  private:
663  friend class QuadratureRuleFactory<ct,d>;
664  PrismQuadratureRule(int p) : QuadratureRule<ct,3>(GeometryType(GeometryType::prism, d))
665  {
666  int m=6;
667  this->delivered_order = PrismQuadraturePointsSingleton<3>::prqp.order(m);
668  for(int i=0; i<m; ++i)
669  {
670  FieldVector<ct,3> local;
671  for (int k=0; k<d; k++)
672  local[k] = PrismQuadraturePointsSingleton<3>::prqp.point(m,i)[k];
673  double weight =
675  // put in container
676  this->push_back(QuadraturePoint<ct,d>(local,weight));
677  }
678  }
679  };
680 
687  template<typename ctype, int dim>
689  private:
690  friend class QuadratureRules<ctype, dim>;
691  static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
692  {
693  return TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
694  }
695  static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
696  {
697  return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
698  }
699  };
700 
701  template<typename ctype>
703  private:
704  enum { dim = 0 };
705  friend class QuadratureRules<ctype, dim>;
706  static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
707  {
708  if (t.isVertex())
709  {
710  return std::numeric_limits<int>::max();
711  }
712  DUNE_THROW(Exception, "Unknown GeometryType");
713  }
714  static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
715  {
716  if (t.isVertex())
717  {
718  return PointQuadratureRule<ctype>();
719  }
720  DUNE_THROW(Exception, "Unknown GeometryType");
721  }
722  };
723 
724  template<typename ctype>
725  class QuadratureRuleFactory<ctype, 1> {
726  private:
727  enum { dim = 1 };
728  friend class QuadratureRules<ctype, dim>;
729  static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
730  {
731  if (t.isLine())
732  {
733  switch (qt) {
734  case QuadratureType::GaussLegendre :
735  return GaussQuadratureRule1D<ctype>::highest_order;
736  case QuadratureType::GaussJacobi_1_0 :
737  return Jacobi1QuadratureRule1D<ctype>::highest_order;
738  case QuadratureType::GaussJacobi_2_0 :
739  return Jacobi2QuadratureRule1D<ctype>::highest_order;
740  case QuadratureType::GaussLobatto :
741  return GaussLobattoQuadratureRule1D<ctype>::highest_order;
742  default :
743  DUNE_THROW(Exception, "Unknown QuadratureType");
744  }
745  }
746  DUNE_THROW(Exception, "Unknown GeometryType");
747  }
748  static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
749  {
750  if (t.isLine())
751  {
752  switch (qt) {
753  case QuadratureType::GaussLegendre :
754  return GaussQuadratureRule1D<ctype>(p);
755  case QuadratureType::GaussJacobi_1_0 :
756  return Jacobi1QuadratureRule1D<ctype>(p);
757  case QuadratureType::GaussJacobi_2_0 :
758  return Jacobi2QuadratureRule1D<ctype>(p);
759  case QuadratureType::GaussLobatto :
760  return GaussLobattoQuadratureRule1D<ctype>(p);
761  default :
762  DUNE_THROW(Exception, "Unknown QuadratureType");
763  }
764  }
765  DUNE_THROW(Exception, "Unknown GeometryType");
766  }
767  };
768 
769  template<typename ctype>
770  class QuadratureRuleFactory<ctype, 2> {
771  private:
772  enum { dim = 2 };
773  friend class QuadratureRules<ctype, dim>;
774  static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
775  {
776  unsigned order =
777  TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
778  if (t.isSimplex())
779  order = std::max
780  (order, unsigned(SimplexQuadratureRule<ctype,dim>::highest_order));
781  return order;
782  }
783  static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
784  {
785  if (t.isSimplex()
786  && qt == QuadratureType::GaussLegendre
787  && p <= SimplexQuadratureRule<ctype,dim>::highest_order)
788  {
789  return SimplexQuadratureRule<ctype,dim>(p);
790  }
791  return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
792  }
793  };
794 
795  template<typename ctype>
796  class QuadratureRuleFactory<ctype, 3> {
797  private:
798  enum { dim = 3 };
799  friend class QuadratureRules<ctype, dim>;
800  static unsigned maxOrder(const GeometryType &t, QuadratureType::Enum qt)
801  {
802  unsigned order =
803  TensorProductQuadratureRule<ctype,dim>::maxOrder(t.id(), qt);
804  if (t.isSimplex())
805  order = std::max
806  (order, unsigned(SimplexQuadratureRule<ctype,dim>::highest_order));
807  if (t.isPrism())
808  order = std::max
809  (order, unsigned(PrismQuadratureRule<ctype,dim>::highest_order));
810  return order;
811  }
812  static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
813  {
814  if (t.isSimplex()
815  && qt == QuadratureType::GaussLegendre
816  && p <= SimplexQuadratureRule<ctype,dim>::highest_order)
817  {
818  return SimplexQuadratureRule<ctype,dim>(p);
819  }
820  if (t.isPrism()
821  && qt == QuadratureType::GaussLegendre
822  && p <= PrismQuadratureRule<ctype,dim>::highest_order)
823  {
824  return PrismQuadratureRule<ctype,dim>(p);
825  }
826  return TensorProductQuadratureRule<ctype,dim>(t.id(), p, qt);
827  }
828  };
829 
830 } // end namespace
831 
832 #endif // DUNE_GEOMETRY_QUADRATURERULES_HH
Jacobi-Gauss quadrature for alpha=2, beta=0.
Definition: quadraturerules.hh:458
Gauss quadrature rule in 1D.
Definition: quadraturerules.hh:265
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
BasicType
Each entity can be tagged by one of these basic types plus its space dimension.
Definition: type.hh:29
@ cube
Cube element in any nonnegative dimension.
Definition: type.hh:31
Jacobi-Gauss quadrature for alpha=1, beta=0.
Definition: quadraturerules.hh:325
Jacobi-Gauss quadrature for alpha=2, beta=0.
Definition: quadraturerules.hh:391
static DUNE_CONSTEXPR std::size_t size(std::size_t dim)
Compute total number of geometry types for the given dimension.
Definition: typeindex.hh:58
static std::size_t index(const GeometryType &gt)
Compute the index for the given geometry type within its dimension.
Definition: typeindex.hh:70
Default exception for dummy implementations.
Definition: exceptions.hh:288
Definition: quadraturerules.hh:515
double weight(int m, int i)
Definition: quadraturerules.hh:606
int order(int m)
Definition: quadraturerules.hh:612
PrismQuadraturePoints()
initialize quadrature points on the interval for all orders
Definition: quadraturerules.hh:521
FieldVector< double, 3 > point(int m, int i)
Definition: quadraturerules.hh:600
Definition: quadraturerules.hh:510
Quadrature rules for prisms.
Definition: quadraturerules.hh:652
Quadrature rules for prisms.
Definition: quadraturerules.hh:645
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
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 Vector & position() const
return local coordinates of integration point i
Definition: quadraturerules.hh:61
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:688
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 &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:225
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:230
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:243
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:132
Dune namespace.
Definition: alignment.hh:10
Standard Dune debug streams.
Definition: quadraturerules.hh:436
Definition: quadraturerules.hh:245
Definition: quadraturerules.hh:303
Definition: quadraturerules.hh:369
Singleton holding the Prism Quadrature points.
Definition: quadraturerules.hh:637
Singleton holding the Prism Quadrature points.
Definition: quadraturerules.hh:629
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.80.0 (May 8, 22:30, 2024)