SecondOrderDivTestvecDivTrialvec.hpp 7.09 KB
Newer Older
1
2
3
4
#pragma once

#include <type_traits>

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

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

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


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

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

  public:
34
35
    GridFunctionOperator(tag::divtestvec_divtrialvec, GridFct const& expr)
      : Super(expr, 2)
36
37
    {}

38
39
    template <class CG, class RN, class CN, class Mat>
    void getElementMatrix(CG const& contextGeo, RN const& rowNode, CN const& colNode, Mat& elementMatrix)
40
    {
41
      static_assert(RN::isPower && CN::isPower,
42
        "Operator can be applied to Power-Nodes only.");
43

44
      const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
45
46
      const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();

47
      auto const& quad = this->getQuadratureRule(contextGeo.type(), rowNode, colNode);
48
      if (sameFE && sameNode)
49
        getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, elementMatrix);
50
      else
51
        getElementMatrixStandard(contextGeo, quad, rowNode, colNode, elementMatrix);
52
53
54
    }

  protected:
55
56
57
58
    template <class CG, class QR, class RN, class CN, class Mat>
    void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
                                  RN const& rowNode, CN const& colNode,
                                  Mat& elementMatrix)
59
    {
60
61
      static const std::size_t CHILDREN = RN::CHILDREN;
      static_assert( RN::CHILDREN == CN::CHILDREN, "" );
62
63
64

      auto const& rowLocalFE = rowNode.child(0).finiteElement();
      auto const& colLocalFE = colNode.child(0).finiteElement();
65
66
      std::size_t rowSize = rowLocalFE.size();
      std::size_t colSize = colLocalFE.size();
67

68
69
      using RowFieldType = typename NodeQuadCache<typename RN::ChildType>::Traits::RangeFieldType;
      using ColFieldType = typename NodeQuadCache<typename CN::ChildType>::Traits::RangeFieldType;
70

71
72
      NodeQuadCache<typename RN::ChildType> rowCache(rowLocalFE.localBasis());
      NodeQuadCache<typename CN::ChildType> colCache(colLocalFE.localBasis());
73

74
75
      auto const& rowGradientsCache = rowCache.evaluateJacobianAtQP(contextGeo, quad);
      auto const& colGradientsCache = colCache.evaluateJacobianAtQP(contextGeo, quad);
76
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
77
        // Position of the current quadrature point in the reference element
78
        decltype(auto) local = contextGeo.local(quad[iq].position());
79
80

        // The transposed inverse Jacobian of the map from the reference element to the element
81
        const auto jacobian = contextGeo.geometry().jacobianInverseTransposed(local);
82
83

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

        // The gradients of the shape functions on the reference element
87
88
        auto const& rowShapeGradients = rowGradientsCache[iq];
        auto const& colShapeGradients = colGradientsCache[iq];
89

90
        // Compute the shape function gradients on the real element
91
        using RowWorldVector = FieldVector<RowFieldType,CG::dow>;
92
        thread_local std::vector<RowWorldVector> rowGradients;
93
        rowGradients.resize(rowShapeGradients.size());
94

95
        for (std::size_t i = 0; i < rowGradients.size(); ++i)
96
          jacobian.mv(rowShapeGradients[i][0], rowGradients[i]);
97

98
        using ColWorldVector = FieldVector<ColFieldType,CG::dow>;
99
        thread_local std::vector<ColWorldVector> colGradients;
100
        colGradients.resize(colShapeGradients.size());
101

102
        for (std::size_t i = 0; i < colGradients.size(); ++i)
103
          jacobian.mv(colShapeGradients[i][0], colGradients[i]);
104

105
106
        for (std::size_t i = 0; i < rowSize; ++i) {
          for (std::size_t j = 0; j < colSize; ++j) {
107
            for (std::size_t k = 0; k < CHILDREN; ++k) {
108
109
              const auto local_ki = rowNode.child(k).localIndex(i);
              const auto local_kj = colNode.child(k).localIndex(j);
110
111
112
113
114
115
116
117
              elementMatrix[local_ki][local_kj] += factor * rowGradients[i][k] * colGradients[j][k];
            }
          }
        }
      }
    }


118
119
120
121
    template <class CG, class QR, class Node, class Mat>
    void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
                                   Node const& node, Node const& /*colNode*/,
                                   Mat& elementMatrix)
122
    {
123
      static const std::size_t CHILDREN = Node::CHILDREN;
124
125

      auto const& localFE = node.child(0).finiteElement();
126
      std::size_t size = localFE.size();
127

128
129
      using LeafNode = typename Node::ChildType;
      using RangeFieldType = typename NodeQuadCache<LeafNode>::Traits::RangeFieldType;
130

131
      NodeQuadCache<LeafNode> cache(localFE.localBasis());
132

133
      auto const& shapeGradientsCache = cache.evaluateJacobianAtQP(contextGeo, quad);
134
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
135
        // Position of the current quadrature point in the reference element
136
        decltype(auto) local = contextGeo.local(quad[iq].position());
137
138

        // The transposed inverse Jacobian of the map from the reference element to the element
139
        const auto jacobian = contextGeo.geometry().jacobianInverseTransposed(local);
140
141

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

        // The gradients of the shape functions on the reference element
145
        auto const& shapeGradients = shapeGradientsCache[iq];
146
147

        // Compute the shape function gradients on the real element
148
        using WorldVector = FieldVector<RangeFieldType,CG::dow>;
149
        thread_local std::vector<WorldVector> gradients;
150
        gradients.resize(shapeGradients.size());
151

152
        for (std::size_t i = 0; i < gradients.size(); ++i)
153
          jacobian.mv(shapeGradients[i][0], gradients[i]);
154
155

        for (std::size_t k = 0; k < CHILDREN; ++k) {
156
          for (std::size_t i = 0; i < size; ++i) {
157
158
            const auto local_ki = node.child(k).localIndex(i);
            const auto value = factor * gradients[i][k];
159
160
            elementMatrix[local_ki][local_ki] += value * gradients[i][k];

161
            for (std::size_t j = i+1; j < size; ++j) {
162
              const auto local_kj = node.child(k).localIndex(j);
163
164
165
166
167
168
169
170
171
              elementMatrix[local_ki][local_kj] += value * gradients[j][k];
              elementMatrix[local_kj][local_ki] += value * gradients[j][k];
            }
          }
        }
      }
    }
  };

172
173
  /** @} **/

174
} // end namespace AMDiS