SecondOrderGradTestGradTrial.hpp 9.79 KB
Newer Older
1
2
3
4
#pragma once

#include <type_traits>

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

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

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


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

32
    using expr_value_type = typename GridFct::Range;
33
34
35
36
    static_assert( Category::Scalar<expr_value_type> || Category::Matrix<expr_value_type>,
      "Expression must be of scalar or matrix type." );

  public:
37
38
    GridFunctionOperator(tag::gradtest_gradtrial, GridFct const& expr)
      : Super(expr, 2)
39
40
    {}

41
42
43
44
45
    template <class Context, class RowNode, class ColNode, class ElementMatrix>
    void getElementMatrix(Context const& context,
                          RowNode const& rowNode,
                          ColNode const& colNode,
                          ElementMatrix& elementMatrix)
46
    {
47
48
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only.");
49

50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
      const bool sameFE = std::is_same<FiniteElementType_t<RowNode>, FiniteElementType_t<ColNode>>::value;
      const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
      using Category = ValueCategory_t<typename GridFct::Range>;

      auto const& quad = this->getQuadratureRule(context.type(), rowNode, colNode);
      if (sameFE && sameNode)
        getElementMatrixOptimized(context, quad, rowNode, colNode, elementMatrix, Category{});
      else if (sameFE)
        getElementMatrixStandard(context, quad, rowNode, colNode, elementMatrix);
      else
        error_exit("Not implemented: currently only the implementation for equal fespaces available");
    }

  protected:

    template <class Context, class QuadratureRule, class RowNode, class ColNode, class ElementMatrix>
    void getElementMatrixStandard(Context const& context,
                                  QuadratureRule const& quad,
                                  RowNode const& rowNode,
                                  ColNode const& colNode,
                                  ElementMatrix& elementMatrix)
    {
      auto const& localFE = rowNode.finiteElement();
73

74
75
      using LocalBasisType = typename std::decay_t<decltype(localFE)>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;
76

77
78
79
80
      LocalBasisCache<LocalBasisType> cache(localFE.localBasis());

      auto const& shapeGradientsCache = cache.evaluateJacobianAtQp(context, quad);
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
81
        // Position of the current quadrature point in the reference element
82
        decltype(auto) local = context.local(quad[iq].position());
83
84

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

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

        // The gradients of the shape functions on the reference element
92
        auto const& shapeGradients = shapeGradientsCache[iq];
93
94

        // Compute the shape function gradients on the real element
95
96
        using WorldVector = Dune::FieldVector<RangeFieldType,Context::dow>;
        thread_local std::vector<WorldVector> gradients;
97
        gradients.resize(shapeGradients.size());
98

99
        for (std::size_t i = 0; i < gradients.size(); ++i)
100
          jacobian.mv(shapeGradients[i][0], gradients[i]);
101

102
103
104
105
        for (std::size_t i = 0; i < localFE.size(); ++i) {
          const auto local_i = rowNode.localIndex(i);
          for (std::size_t j = 0; j < localFE.size(); ++j) {
            const auto local_j = colNode.localIndex(j);
106
107
108
109
110
111
            elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
          }
        }
      }
    }

112
113
114
115
116
117
118
    template <class Context, class QuadratureRule, class Node, class ElementMatrix>
    void getElementMatrixOptimized(Context const& context,
                                   QuadratureRule const& quad,
                                   Node const& node,
                                   Node const& /*colNode*/,
                                   ElementMatrix& elementMatrix,
                                   tag::scalar)
119
120
121
    {
      auto const& localFE = node.finiteElement();

122
123
124
      using LocalBasisType = typename std::decay_t<decltype(localFE)>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;

125
      LocalBasisCache<LocalBasisType> cache(localFE.localBasis());
126
      auto const& shapeGradientsCache = cache.evaluateJacobianAtQp(context, quad);
127
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
128
        // Position of the current quadrature point in the reference element
129
        decltype(auto) local = context.local(quad[iq].position());
130
131

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

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

        // The gradients of the shape functions on the reference element
139
        auto const& shapeGradients = shapeGradientsCache[iq];
140
141

        // Compute the shape function gradients on the real element
142
143
        using WorldVector = Dune::FieldVector<RangeFieldType,Context::dow>;
        thread_local std::vector<WorldVector> gradients;
144
        gradients.resize(shapeGradients.size());
145
        for (std::size_t i = 0; i < gradients.size(); ++i)
146
          jacobian.mv(shapeGradients[i][0], gradients[i]);
147
148

        for (std::size_t i = 0; i < localFE.size(); ++i) {
149
          const auto local_i = node.localIndex(i);
150
151
152
153

          elementMatrix[local_i][local_i] += factor * (gradients[i] * gradients[i]);

          for (std::size_t j = i+1; j < localFE.size(); ++j) {
154
155
            const auto local_j = node.localIndex(j);
            const auto value = factor * (gradients[i] * gradients[j]);
156
157
158
159
160
161
162
163

            elementMatrix[local_i][local_j] += value;
            elementMatrix[local_j][local_i] += value;
          }
        }
      }
    }

164
165
166
167
168
169
170
    template <class Context, class QuadratureRule, class Node, class ElementMatrix>
    void getElementMatrixOptimized(Context const& context,
                                   QuadratureRule const& quad,
                                   Node const& node,
                                   Node const& /*colNode*/,
                                   ElementMatrix& elementMatrix,
                                   tag::matrix)
171
172
173
    {
      auto const& localFE = node.finiteElement();

174
175
176
      using LocalBasisType = typename std::decay_t<decltype(localFE)>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;

177
178
179
      LocalBasisCache<LocalBasisType> cache(localFE.localBasis());
      auto const& shapeGradientsCache = cache.evaluateJacobianAtQp(context, quad);
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
180
        // Position of the current quadrature point in the reference element
181
        decltype(auto) local = context.local(quad[iq].position());
182
183

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

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

        // The gradients of the shape functions on the reference element
191
        auto const& shapeGradients = shapeGradientsCache[iq];
192
193

        // Compute the shape function gradients on the real element
194
195
        using WorldVector = Dune::FieldVector<RangeFieldType,Context::dow>;
        thread_local std::vector<WorldVector> gradients;
196
        gradients.resize(shapeGradients.size());
197
        for (std::size_t i = 0; i < gradients.size(); ++i)
198
          jacobian.mv(shapeGradients[i][0], gradients[i]);
199
200

        for (std::size_t i = 0; i < localFE.size(); ++i) {
201
          const auto local_i = node.localIndex(i);
202
          for (std::size_t j = 0; j < localFE.size(); ++j) {
203
            const auto local_j = node.localIndex(j);
204
205
206
207
208
209
210
211
            elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
          }
        }
      }
    }

  protected:

212
213
214
215
216
    template <class S, class F, class T, int dow,
      std::enable_if_t<Category::Scalar<S>,int> = 0>
    T eval(S const& scalar, F factor,
           Dune::FieldVector<T,dow> const& grad_test,
           Dune::FieldVector<T,dow> const& grad_trial) const
217
218
219
220
    {
      return (scalar * factor) * (grad_test * grad_trial);
    }

221
222
223
224
225
    template <class M, class F, class T, int dow,
      std::enable_if_t<Category::Matrix<M>,int> = 0>
    T eval(M const& mat, F factor,
           Dune::FieldVector<T,dow> const& grad_test,
           Dune::FieldVector<T,dow> const& grad_trial) const
226
    {
227
      return factor * (grad_test * (mat * grad_trial));
228
229
230
    }
  };

231
232
  /** @} **/

233
} // end namespace AMDiS