AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
QuadratureFactory.hpp
1 #pragma once
2 
3 #include <dune/common/exceptions.hh>
4 #include <dune/geometry/quadraturerules.hh>
5 #include <dune/geometry/type.hh>
6 
7 #include <amdis/common/Order.hpp>
8 #include <amdis/functions/Order.hpp>
9 
10 namespace AMDiS
11 {
12  template <class Geometry, class RN, class CN>
13  int getQuadratureDegree(Geometry const& geometry, int derivOrder, int coeffDegree,
14  RN const& rowNode, CN const& colNode)
15  {
16  if (coeffDegree < 0)
17  DUNE_THROW(Dune::Exception, "Polynomial degree of coefficients cannot be determined. "
18  "Please provide a quadrature order manually.");
19 
20  int psiDegree = order(rowNode);
21  int phiDegree = order(colNode);
22 
23  int degree = psiDegree + phiDegree + coeffDegree;
24  if (geometry.type().isSimplex())
25  degree -= derivOrder;
26  if (geometry.affine())
27  degree += 1;
28 
29  return degree;
30  }
31 
32  template <class Geometry, class Node>
33  int getQuadratureDegree(Geometry const& geometry, int derivOrder, int coeffDegree,
34  Node const& node)
35  {
36  if (coeffDegree < 0)
37  DUNE_THROW(Dune::Exception, "Polynomial degree of GridFunction cannot be determined. "
38  "Please provide a quadrature order manually.");
39 
40  int psiDegree = order(node);
41  int degree = psiDegree + coeffDegree;
42 
43  if (geometry.type().isSimplex())
44  degree -= derivOrder;
45  if (geometry.affine())
46  degree += 1;
47 
48  return degree;
49  }
50 
51  template <class Geometry, class... Nodes>
52  auto const& getQuadratureRule(Geometry const& geometry, int derivOrder,
53  int coeffDegree, Nodes const&... nodes)
54  {
55  int degree = getQuadratureDegree(geometry, derivOrder, coeffDegree, nodes...);
56  using QuadRules
57  = Dune::QuadratureRules<typename Geometry::ctype, Geometry::mydimension>;
58  return QuadRules::rule(geometry.type(), degree);
59  }
60 
61 } // end namespace AMDiS
Definition: AdaptBase.hpp:6