AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
ZeroOrderTestvec.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 testvec {};
18  }
19 
20 
23  {
24  public:
26 
27  template <class CG, class Node, class Quad, class LocalFct, class Vec>
28  void assemble(CG const& contextGeo, Node const& node, Quad const& quad,
29  LocalFct const& localFct, Vec& elementVector) const
30  {
31  static_assert(static_size_v<typename LocalFct::Range> == CG::dow,
32  "Expression must be of vector type." );
33  static_assert(Node::isPower,
34  "Operator can be applied to Power-Nodes only.");
35 
36  static const std::size_t CHILDREN = Node::CHILDREN;
37  static_assert(CHILDREN == CG::dow);
38 
39  std::size_t size = node.child(0).size();
40 
41  for (auto const& qp : quad) {
42  // Position of the current quadrature point in the reference element
43  auto&& local = contextGeo.local(qp.position());
44 
45  // The multiplicative factor in the integral transformation formula
46  const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
47  const auto exprValue = localFct(local);
48 
49  auto const& shapeValues = node.child(0).localBasisValuesAt(local);
50 
51  for (std::size_t i = 0; i < size; ++i) {
52  const auto value = exprValue * (factor * shapeValues[i]);
53  for (std::size_t k = 0; k < CHILDREN; ++k) {
54  const auto local_ki = node.child(k).localIndex(i);
55  elementVector[local_ki] += at(value,k);
56  }
57  }
58  }
59  }
60  };
61 
62  template <class LC>
63  struct GridFunctionOperatorRegistry<tag::testvec, LC>
64  {
65  static constexpr int degree = 0;
66  using type = ZeroOrderTestVec;
67  };
68 
71 } // end namespace AMDiS
Definition: AdaptBase.hpp:6
zero-order vector-operator
Definition: ZeroOrderTestvec.hpp:22
Definition: ZeroOrderTestvec.hpp:17
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:204