AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
FirstOrderDivTestvec.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 divtestvec {};
18  }
19 
20 
23  {
24  public:
26 
27  template <class CG, class Node, class Quad, class LF, class Vec>
28  void assemble(CG const& context, Node const& node, Quad const& quad,
29  LF const& localFct, Vec& elementVector) const
30  {
31  static_assert(Node::isPower, "Node must be a Power-Node.");
32  static_assert(static_size_v<typename LF::Range> == 1,
33  "Expression must be of scalar type.");
34  static_assert(Node::CHILDREN == CG::dow);
35 
36  std::size_t feSize = node.child(0).size();
37 
38  using RangeFieldType = typename Node::ChildType::LocalBasis::Traits::RangeFieldType;
39  using WorldVector = FieldVector<RangeFieldType,CG::dow>;
40  std::vector<WorldVector> gradients;
41 
42  for (auto const& qp : quad) {
43  // Position of the current quadrature point in the reference element
44  auto&& local = context.local(qp.position());
45 
46  // The transposed inverse Jacobian of the map from the reference element to the element
47  const auto jacobian = context.geometry().jacobianInverseTransposed(local);
48 
49  // The multiplicative factor in the integral transformation formula
50  const auto factor = localFct(local) * context.integrationElement(qp.position()) * qp.weight();
51 
52  // The gradients of the shape functions on the reference element
53  auto const& shapeGradients = node.child(0).localBasisJacobiansAt(local);
54 
55  // Compute the shape function gradients on the real element
56  gradients.resize(shapeGradients.size());
57 
58  for (std::size_t i = 0; i < gradients.size(); ++i)
59  jacobian.mv(shapeGradients[i][0], gradients[i]);
60 
61  for (std::size_t j = 0; j < feSize; ++j) {
62  for (std::size_t k = 0; k < CG::dow; ++k) {
63  const auto local_kj = node.child(k).localIndex(j);
64  elementVector[local_kj] += factor * gradients[j][k];
65  }
66  }
67  }
68 
69  }
70  };
71 
72  template <class LC>
73  struct GridFunctionOperatorRegistry<tag::divtestvec, LC>
74  {
75  static constexpr int degree = 1;
76 
77  using type = FirstOrderDivTestVec;
78  };
79 
82 } // end namespace AMDiS
Definition: FirstOrderDivTestvec.hpp:17
Definition: AdaptBase.hpp:6
first-order operator
Definition: FirstOrderDivTestvec.hpp:22
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:204