AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
FirstOrderTestDivTrialvec.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_divtrialvec {};
18  }
19 
20 
23  {
24  public:
26 
27  template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
28  void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
29  Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
30  {
31  static_assert(static_size_v<typename LocalFct::Range> == 1, "Expression must be of scalar type.");
32  static_assert(RN::isLeaf && CN::isPower,
33  "row-node must be Leaf-Node and col-node must be a Power-Node.");
34 
35  assert( colNode.degree() == CG::dow );
36 
37  std::size_t rowSize = rowNode.size();
38  std::size_t colSize = colNode.child(0).size();
39 
40  using RangeFieldType = typename CN::ChildType::LocalBasis::Traits::RangeFieldType;
41  using WorldVector = FieldVector<RangeFieldType,CG::dow>;
42  std::vector<WorldVector> colGradients;
43 
44  for (auto const& qp : quad) {
45  // Position of the current quadrature point in the reference element
46  auto&& local = contextGeo.coordinateInElement(qp.position());
47 
48  // The transposed inverse Jacobian of the map from the reference element to the element
49  const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
50 
51  // The multiplicative factor in the integral transformation formula
52  const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
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.child(0).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 i = 0; i < rowSize; ++i) {
67  const auto local_i = rowNode.localIndex(i);
68  const auto value = factor * shapeValues[i];
69  for (std::size_t j = 0; j < colSize; ++j) {
70  for (std::size_t k = 0; k < colNode.degree(); ++k) {
71  const auto local_kj = colNode.child(k).localIndex(j);
72  elementMatrix[local_i][local_kj] += value * colGradients[j][k];
73  }
74  }
75  }
76  }
77 
78  }
79  };
80 
81  template <class LC>
82  struct GridFunctionOperatorRegistry<tag::test_divtrialvec, LC>
83  {
84  static constexpr int degree = 1;
86  };
87 
90 } // end namespace AMDiS
Definition: AdaptBase.hpp:6
first-order operator
Definition: FirstOrderTestDivTrialvec.hpp:22
Definition: FirstOrderTestDivTrialvec.hpp:17
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:216