AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
FirstOrderTestPartialTrial.hpp
1 #pragma once
2 
3 #include <type_traits>
4 
5 #include <amdis/GridFunctionOperator.hpp>
6 #include <amdis/Output.hpp>
7 #include <amdis/common/StaticSize.hpp>
8 #include <amdis/utility/LocalToGlobalAdapter.hpp>
9 
10 namespace AMDiS
11 {
17  namespace tag
18  {
20  {
21  std::size_t comp;
22  };
23 
25  {
26  std::size_t comp;
27  };
28  }
29 
30 
32  template <class LC, class GridFct>
33  class GridFunctionOperator<tag::test_partialtrial, LC, GridFct>
34  : public GridFunctionOperatorBase<GridFunctionOperator<tag::test_partialtrial, LC, GridFct>, LC, GridFct>
35  {
36  using Self = GridFunctionOperator;
38 
39  static_assert( static_size_v<typename GridFct::Range> == 1, "Expression must be of scalar type." );
40 
41  public:
43  : Super(expr, 1)
44  , comp_(tag.comp)
45  {}
46 
47  template <class CG, class RN, class CN, class Mat>
48  void getElementMatrix(CG const& contextGeo, RN const& rowNode, CN const& colNode, Mat& elementMatrix)
49  {
50  static_assert(RN::isLeaf && CN::isLeaf,
51  "Operator can be applied to Leaf-Nodes only.");
52 
53  auto const& quad = this->getQuadratureRule(contextGeo.type(), rowNode, colNode);
54 
55  std::size_t rowSize = rowNode.size();
56  std::size_t colSize = colNode.size();
57 
58  LocalToGlobalBasisAdapter colLocalBasis(colNode, contextGeo.geometry());
59  for (auto const& qp : quad) {
60  // Position of the current quadrature point in the reference element
61  decltype(auto) local = contextGeo.local(qp.position());
62 
63  // The multiplicative factor in the integral transformation formula
64  auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
65  factor *= Super::coefficient(local);
66 
67  // the values of the shape functions on the reference element at the quadrature point
68  auto const& shapeValues = rowNode.localBasisValuesAt(local);
69 
70  // Compute the shape function gradients on the real element
71  auto const& colPartial = colLocalBasis.partialsAt(local, comp_);
72 
73  for (std::size_t j = 0; j < colSize; ++j) {
74  const auto local_j = colNode.localIndex(j);
75  const auto value = factor * colPartial[j];
76  for (std::size_t i = 0; i < rowSize; ++i) {
77  const auto local_i = rowNode.localIndex(i);
78  elementMatrix[local_i][local_j] += value * shapeValues[i];
79  }
80  }
81  }
82  }
83 
84  private:
85  std::size_t comp_;
86  };
87 
88 
90  template <class Expr, class... QuadratureArgs>
91  auto fot(Expr&& expr, tag::partial_trial t, QuadratureArgs&&... args)
92  {
93  return makeOperator(tag::test_partialtrial{t.comp}, FWD(expr), FWD(args)...);
94  }
95 
98 } // 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
Convert a simple (scalar) local basis into a global basis.
Definition: LocalToGlobalAdapter.hpp:64
Definition: FirstOrderTestPartialTrial.hpp:24
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: FirstOrderTestPartialTrial.hpp:19
auto fot(Expr &&expr, tag::grad_test, QuadratureArgs &&... args)
Create a first-order term with derivative on trial-function.
Definition: FirstOrderGradTestTrial.hpp:40
The main implementation of an CRTP-base class for operators using a grid-function coefficient to be u...
Definition: GridFunctionOperator.hpp:39