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/Output.hpp>
7 #include <amdis/common/StaticSize.hpp>
8 
9 namespace AMDiS
10 {
16  namespace tag
17  {
18  struct divtestvec {};
19  }
20 
21 
23  template <class LC, class GridFct>
24  class GridFunctionOperator<tag::divtestvec, LC, GridFct>
25  : public GridFunctionOperatorBase<GridFunctionOperator<tag::divtestvec, LC, GridFct>, LC, GridFct>
26  {
27  using Self = GridFunctionOperator;
29 
30  static_assert( static_size_v<typename GridFct::Range> == 1, "Expression must be of scalar type." );
31 
32  public:
34  : Super(expr, 1)
35  {}
36 
37  template <class CG, class Node, class Vec>
38  void getElementVector(CG const& context, Node const& node, Vec& elementVector)
39  {
40  static_assert(Node::isPower, "Node must be a Power-Node.");
41 
42  static const std::size_t CHILDREN = Node::CHILDREN;
43  static_assert( CHILDREN == CG::dow, "" );
44 
45  auto const& quad = this->getQuadratureRule(context.type(), node);
46  std::size_t feSize = node.child(0).size();
47 
48  using RangeFieldType = typename Node::ChildType::LocalBasis::Traits::RangeFieldType;
49  using WorldVector = FieldVector<RangeFieldType,CG::dow>;
50  std::vector<WorldVector> gradients;
51 
52  for (auto const& qp : quad) {
53  // Position of the current quadrature point in the reference element
54  auto&& local = context.local(qp.position());
55 
56  // The transposed inverse Jacobian of the map from the reference element to the element
57  const auto jacobian = context.geometry().jacobianInverseTransposed(local);
58 
59  // The multiplicative factor in the integral transformation formula
60  const auto factor = Super::coefficient(local) * context.integrationElement(qp.position()) * qp.weight();
61 
62  // The gradients of the shape functions on the reference element
63  auto const& shapeGradients = node.child(0).localBasisJacobiansAt(local);
64 
65  // Compute the shape function gradients on the real element
66  gradients.resize(shapeGradients.size());
67 
68  for (std::size_t i = 0; i < gradients.size(); ++i)
69  jacobian.mv(shapeGradients[i][0], gradients[i]);
70 
71  for (std::size_t j = 0; j < feSize; ++j) {
72  for (std::size_t k = 0; k < CHILDREN; ++k) {
73  const auto local_kj = node.child(k).localIndex(j);
74  elementVector[local_kj] += factor * gradients[j][k];
75  }
76  }
77  }
78 
79  }
80  };
81 
84 } // end namespace AMDiS
The base-template for GridFunctionOperators.
Definition: GridFunctionOperator.hpp:242
Definition: FirstOrderDivTestvec.hpp:18
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
The main implementation of an CRTP-base class for operators using a grid-function coefficient to be u...
Definition: GridFunctionOperator.hpp:39