AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
FirstOrderTestGradTrial.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_gradtrial {};
18  struct grad_trial {};
19  }
20 
21 
24  {
25  public:
27 
28  template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
29  void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
30  Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
31  {
32  static_assert(static_size_v<typename LocalFct::Range> == CG::dow, "Expression must be of vector type.");
33  static_assert(RN::isLeaf && CN::isLeaf,
34  "Operator can be applied to Leaf-Nodes only.");
35 
36  std::size_t rowSize = rowNode.size();
37  std::size_t colSize = colNode.size();
38 
39  using RangeFieldType = typename CN::LocalBasis::Traits::RangeFieldType;
40  using WorldVector = FieldVector<RangeFieldType,CG::dow>;
41  std::vector<WorldVector> colGradients;
42 
43  for (auto const& qp : quad) {
44  // Position of the current quadrature point in the reference element
45  auto&& local = contextGeo.coordinateInElement(qp.position());
46 
47  // The transposed inverse Jacobian of the map from the reference element to the element
48  const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
49 
50  // The multiplicative factor in the integral transformation formula
51  const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
52  const auto b = localFct(local);
53 
54  // the values of the shape functions on the reference element at the quadrature point
55  auto const& shapeValues = rowNode.localBasisValuesAt(local);
56 
57  // The gradients of the shape functions on the reference element
58  auto const& shapeGradients = colNode.localBasisJacobiansAt(local);
59 
60  // Compute the shape function gradients on the real element
61  colGradients.resize(shapeGradients.size());
62 
63  for (std::size_t i = 0; i < colGradients.size(); ++i)
64  jacobian.mv(shapeGradients[i][0], colGradients[i]);
65 
66  for (std::size_t j = 0; j < colSize; ++j) {
67  const auto local_j = colNode.localIndex(j);
68  const auto value = factor * (b * colGradients[j]);
69  for (std::size_t i = 0; i < rowSize; ++i) {
70  const auto local_i = rowNode.localIndex(i);
71  elementMatrix[local_i][local_j] += value * shapeValues[i];
72  }
73  }
74  }
75 
76  }
77  };
78 
79  template <class LC>
80  struct GridFunctionOperatorRegistry<tag::test_gradtrial, LC>
81  {
82  static constexpr int degree = 1;
84  };
85 
86 
88  template <class Expr>
89  auto fot(Expr&& expr, tag::grad_trial, int quadOrder = -1)
90  {
91  return makeOperator(tag::test_gradtrial{}, FWD(expr), quadOrder);
92  }
93 
96 } // end namespace AMDiS
Definition: FirstOrderTestGradTrial.hpp:17
auto fot(Expr &&expr, tag::grad_trial, int quadOrder=-1)
Create a first-order term with derivative on test-function.
Definition: FirstOrderTestGradTrial.hpp:89
Definition: AdaptBase.hpp:6
first-order operator
Definition: FirstOrderTestGradTrial.hpp:23
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition: GridFunctionOperator.hpp:235
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:216
Definition: FirstOrderTestGradTrial.hpp:18