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

#include <type_traits>

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

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

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


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

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

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

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

49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
      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();
72

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

77
      for (auto const& qp : quad) {
78
        // Position of the current quadrature point in the reference element
79
        decltype(auto) local = context.local(qp.position());
80
81

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

        // The multiplicative factor in the integral transformation formula
85
        const auto factor = context.integrationElement(qp.position()) * qp.weight();
86
        const auto exprValue = Super::coefficient(local);
87
88

        // The gradients of the shape functions on the reference element
89
90
        thread_local std::vector<JacobianType> shapeGradients;
        localFE.localBasis().evaluateJacobian(local, shapeGradients);
91
92

        // Compute the shape function gradients on the real element
93
94
        thread_local std::vector<Dune::FieldVector<RangeFieldType,Context::dow> > gradients;
        gradients.resize(shapeGradients.size());
95
        for (std::size_t i = 0; i < gradients.size(); ++i)
96
          jacobian.mv(shapeGradients[i][0], gradients[i]);
97

98
99
100
101
        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);
102
103
104
105
106
107
            elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
          }
        }
      }
    }

108
109
110
111
112
113
114
    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)
115
116
117
    {
      auto const& localFE = node.finiteElement();

118
119
120
      using LocalBasisType = typename std::decay_t<decltype(localFE)>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;

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

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

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

        // The gradients of the shape functions on the reference element
135
        auto const& shapeGradients = shapeGradientsCache[iq];
136
137

        // Compute the shape function gradients on the real element
138
139
        thread_local std::vector<Dune::FieldVector<RangeFieldType,Context::dow> > gradients;
        gradients.resize(shapeGradients.size());
140
        for (std::size_t i = 0; i < gradients.size(); ++i)
141
          jacobian.mv(shapeGradients[i][0], gradients[i]);
142
143

        for (std::size_t i = 0; i < localFE.size(); ++i) {
144
          const auto local_i = node.localIndex(i);
145
146
147
148

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

          for (std::size_t j = i+1; j < localFE.size(); ++j) {
149
150
            const auto local_j = node.localIndex(j);
            const auto value = factor * (gradients[i] * gradients[j]);
151
152
153
154
155
156
157
158

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

159
160
161
162
163
164
165
    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)
166
167
168
    {
      auto const& localFE = node.finiteElement();

169
170
171
172
173
      using LocalBasisType = typename std::decay_t<decltype(localFE)>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;
      using JacobianType = typename LocalBasisType::Traits::JacobianType;

      for (auto const& qp : quad) {
174
        // Position of the current quadrature point in the reference element
175
        decltype(auto) local = context.local(qp.position());
176
177

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

        // The multiplicative factor in the integral transformation formula
181
        const auto factor = context.integrationElement(qp.position()) * qp.weight();
182
        const auto exprValue = Super::coefficient(local);
183
184

        // The gradients of the shape functions on the reference element
185
186
        thread_local std::vector<JacobianType> shapeGradients;
        localFE.localBasis().evaluateJacobian(local, shapeGradients);
187
188

        // Compute the shape function gradients on the real element
189
190
        thread_local std::vector<Dune::FieldVector<RangeFieldType,Context::dow> > gradients;
        gradients.resize(shapeGradients.size());
191
        for (std::size_t i = 0; i < gradients.size(); ++i)
192
          jacobian.mv(shapeGradients[i][0], gradients[i]);
193
194

        for (std::size_t i = 0; i < localFE.size(); ++i) {
195
          const auto local_i = node.localIndex(i);
196
          for (std::size_t j = 0; j < localFE.size(); ++j) {
197
            const auto local_j = node.localIndex(j);
198
199
200
201
202
203
204
205
            elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
          }
        }
      }
    }

  protected:

206
207
208
209
210
    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
211
212
213
214
    {
      return (scalar * factor) * (grad_test * grad_trial);
    }

215
216
217
218
219
    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
220
    {
221
      return factor * (grad_test * (mat * grad_trial));
222
223
224
    }
  };

225
226
  /** @} **/

227
} // end namespace AMDiS