FirstOrderTestDivTrialvec.hpp 3.92 KB
Newer Older
1
2
3
4
#pragma once

#include <type_traits>

5
#include <amdis/GridFunctionOperator.hpp>
6
#include <amdis/LocalBasisCache.hpp>
7
#include <amdis/Output.hpp>
8
#include <amdis/common/Size.hpp>
9
10
11

namespace AMDiS
{
12
13
14
15
16
  /**
   * \addtogroup operators
   * @{
   **/

17
18
19
20
21
22
  namespace tag
  {
    struct test_divtrialvec {};
  }


23
  /// first-order operator \f$ \langle\psi, c\,\nabla\cdot\Phi\rangle \f$
24
25
26
27
  template <class LocalContext, class GridFct>
  class GridFunctionOperator<tag::test_divtrialvec, LocalContext, GridFct>
      : public GridFunctionOperatorBase<GridFunctionOperator<tag::test_divtrialvec, LocalContext, GridFct>,
                                        LocalContext, GridFct>
28
  {
29
    using Self = GridFunctionOperator;
30
    using Super = GridFunctionOperatorBase<Self, LocalContext, GridFct>;
31

32
    static_assert( Size_v<typename GridFct::Range> == 1, "Expression must be of scalar type." );
33
34

  public:
35
36
    GridFunctionOperator(tag::test_divtrialvec, GridFct const& expr)
      : Super(expr, 1)
37
38
    {}

39
40
41
42
43
    template <class Context, class RowNode, class ColNode, class ElementMatrix>
    void getElementMatrix(Context const& context,
                          RowNode const& rowNode,
                          ColNode const& colNode,
                          ElementMatrix& elementMatrix)
44
    {
45
46
      static_assert(RowNode::isLeaf && ColNode::isPower,
        "RowNode must be Leaf-Node and ColNode must be a Power-Node.");
47
48
49
50

      static const std::size_t CHILDREN = ColNode::CHILDREN;
      static_assert( CHILDREN == Context::dow, "" );

51
      auto const& quad = this->getQuadratureRule(context.type(), rowNode, colNode);
52
53
      auto const& rowLocalFE = rowNode.finiteElement();
      auto const& colLocalFE = colNode.child(0).finiteElement();
54
55
      std::size_t rowSize = rowLocalFE.size();
      std::size_t colSize = colLocalFE.size();
56

57
      using RangeFieldType = typename NodeQuadCache<typename ColNode::ChildType>::Traits::RangeFieldType;
58

59
60
      NodeQuadCache<RowNode> rowCache(rowLocalFE.localBasis());
      NodeQuadCache<typename ColNode::ChildType> colCache(colLocalFE.localBasis());
61

62
63
      auto const& shapeValuesCache = rowCache.evaluateFunctionAtQP(context, quad);
      auto const& shapeGradientsCache = colCache.evaluateJacobianAtQP(context, quad);
64
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
65
        // Position of the current quadrature point in the reference element
66
        decltype(auto) local = context.local(quad[iq].position());
67
68

        // The transposed inverse Jacobian of the map from the reference element to the element
69
        const auto jacobian = context.geometry().jacobianInverseTransposed(local);
70
71

        // The multiplicative factor in the integral transformation formula
72
        const auto factor = Super::coefficient(local) * context.integrationElement(quad[iq].position()) * quad[iq].weight();
73

74
75
        // the values of the shape functions on the reference element at the quadrature point
        auto const& shapeValues = shapeValuesCache[iq];
76
77

        // The gradients of the shape functions on the reference element
78
        auto const& shapeGradients = shapeGradientsCache[iq];
79
80

        // Compute the shape function gradients on the real element
81
82
        using WorldVector = FieldVector<RangeFieldType,Context::dow>;
        thread_local std::vector<WorldVector> colGradients;
83
        colGradients.resize(shapeGradients.size());
84
85

        for (std::size_t i = 0; i < colGradients.size(); ++i)
86
          jacobian.mv(shapeGradients[i][0], colGradients[i]);
87

88
        for (std::size_t i = 0; i < rowSize; ++i) {
89
90
          const auto local_i = rowNode.localIndex(i);
          const auto value = factor * shapeValues[i];
91
          for (std::size_t j = 0; j < colSize; ++j) {
92
            for (std::size_t k = 0; k < CHILDREN; ++k) {
93
              const auto local_kj = colNode.child(k).localIndex(j);
94
95
96
97
98
99
100
101
102
              elementMatrix[local_i][local_kj] += value * colGradients[j][k];
            }
          }
        }
      }

    }
  };

103
104
  /** @} **/

105
} // end namespace AMDiS