AMDiS  0.1
The Adaptive Multi-Dimensional Simulation Toolbox
ZeroOrderTestTrial.hpp
1 #pragma once
2 
3 #include <type_traits>
4 
5 #include <amdis/GridFunctionOperator.hpp>
6 #include <amdis/common/StaticSize.hpp>
7 
8 namespace AMDiS
9 {
15  namespace tag
16  {
17  struct test_trial {};
18  }
19 
20 
22  template <class LC, class GridFct>
23  class GridFunctionOperator<tag::test_trial, LC, GridFct>
24  : public GridFunctionOperatorBase<GridFunctionOperator<tag::test_trial, LC, GridFct>, LC, GridFct>
25  {
26  using Self = GridFunctionOperator;
28 
29  static_assert( static_size_v<typename GridFct::Range> == 1, "Expression must be of scalar type." );
30 
31  public:
33  : Super(expr, 0)
34  {}
35 
36  template <class CG, class RN, class CN, class Mat>
37  void getElementMatrix(CG const& contextGeo, RN const& rowNode, CN const& colNode, Mat& elementMatrix)
38  {
39  const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
40  const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
41 
42  auto const& quad = this->getQuadratureRule(contextGeo.type(), rowNode, colNode);
43  if (sameFE && sameNode)
44  getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, elementMatrix);
45  else
46  getElementMatrixStandard(contextGeo, quad, rowNode, colNode, elementMatrix);
47  }
48 
49  template <class CG, class Node, class Vec>
50  void getElementVector(CG const& contextGeo, Node const& node, Vec& elementVector)
51  {
52  static_assert(Node::isLeaf, "Operator can be applied to Leaf-Nodes only");
53 
54  auto const& quad = this->getQuadratureRule(contextGeo.type(), node);
55  std::size_t size = node.size();
56 
57  for (auto const& qp : quad) {
58  // Position of the current quadrature point in the reference element
59  auto&& local = contextGeo.local(qp.position());
60 
61  // The multiplicative factor in the integral transformation formula
62  const auto factor = Super::coefficient(local) * contextGeo.integrationElement(qp.position())
63  * qp.weight();
64 
65  auto const& shapeValues = node.localBasisValuesAt(local);
66  for (std::size_t i = 0; i < size; ++i) {
67  const auto local_i = node.localIndex(i);
68  elementVector[local_i] += factor * shapeValues[i];
69  }
70  }
71  }
72 
73 
74  protected:
75  template <class CG, class QR, class RN, class CN, class Mat>
76  void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
77  RN const& rowNode, CN const& colNode,
78  Mat& elementMatrix)
79  {
80  static_assert(RN::isLeaf && CN::isLeaf,
81  "Operator can be applied to Leaf-Nodes only.");
82 
83  std::size_t rowSize = rowNode.size();
84  std::size_t colSize = colNode.size();
85 
86  for (auto const& qp : quad) {
87  // Position of the current quadrature point in the reference element
88  auto&& local = contextGeo.local(qp.position());
89 
90  // The multiplicative factor in the integral transformation formula
91  const auto factor = Super::coefficient(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
92 
93  auto const& rowShapeValues = rowNode.localBasisValuesAt(local);
94  auto const& colShapeValues = colNode.localBasisValuesAt(local);
95 
96  for (std::size_t i = 0; i < rowSize; ++i) {
97  const auto local_i = rowNode.localIndex(i);
98  const auto value = factor * rowShapeValues[i];
99 
100  for (std::size_t j = 0; j < colSize; ++j) {
101  const auto local_j = colNode.localIndex(j);
102  elementMatrix[local_i][local_j] += value * colShapeValues[j];
103  }
104  }
105  }
106  }
107 
108 
109  template <class CG, class QR, class RN, class CN, class Mat>
110  void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
111  RN const& node, CN const& /*colNode*/,
112  Mat& elementMatrix)
113  {
114  static_assert(RN::isLeaf && CN::isLeaf,
115  "Operator can be applied to Leaf-Nodes only.");
116 
117  std::size_t size = node.size();
118 
119  for (auto const& qp : quad) {
120  // Position of the current quadrature point in the reference element
121  auto&& local = contextGeo.local(qp.position());
122 
123  // The multiplicative factor in the integral transformation formula
124  const auto factor = Super::coefficient(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
125 
126  auto const& shapeValues = node.localBasisValuesAt(local);
127 
128  for (std::size_t i = 0; i < size; ++i) {
129  const auto local_i = node.localIndex(i);
130 
131  const auto value = factor * shapeValues[i];
132  elementMatrix[local_i][local_i] += value * shapeValues[i];
133 
134  for (std::size_t j = i+1; j < size; ++j) {
135  const auto local_j = node.localIndex(j);
136 
137  elementMatrix[local_i][local_j] += value * shapeValues[j];
138  elementMatrix[local_j][local_i] += value * shapeValues[j];
139  }
140  }
141  }
142  }
143  };
144 
145 
147  template <class Expr, class... QuadratureArgs>
148  auto zot(Expr&& expr, QuadratureArgs&&... args)
149  {
150  return makeOperator(tag::test_trial{}, FWD(expr), FWD(args)...);
151  }
152 
155 } // end namespace AMDiS
auto makeOperator(Tag tag, Expr &&expr, QuadratureArgs &&... args)
Store tag and expression into a PreGridFunctionOperator to create a GridFunctionOperator.
Definition: GridFunctionOperator.hpp:220
The base-template for GridFunctionOperators.
Definition: GridFunctionOperator.hpp:242
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
Definition: ZeroOrderTestTrial.hpp:17
auto zot(Expr &&expr, QuadratureArgs &&... args)
Create a zero-order term.
Definition: ZeroOrderTestTrial.hpp:148
The main implementation of an CRTP-base class for operators using a grid-function coefficient to be u...
Definition: GridFunctionOperator.hpp:39