SecondOrderGradTestGradTrial.hpp 8.62 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
23
24
25
  /// second-order operator \f$ \langle\nabla\psi, c\,\nabla\phi\rangle \f$, or \f$ \langle\nabla\psi, A\,\nabla\phi\rangle \f$
  template <class GridFct, class QuadCreator>
  class GridFunctionOperator<tag::gradtest_gradtrial, GridFct, QuadCreator>
      : public GridFunctionOperatorBase<GridFct, QuadCreator>
26
  {
27
    using Super = GridFunctionOperatorBase<GridFct, QuadCreator>;
28

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

  public:
34
35
    GridFunctionOperator(tag::gradtest_gradtrial, GridFct const& expr, QuadCreator const& quadCreator)
      : Super(expr, quadCreator, 2)
36
37
38
39
40
41
42
43
44
45
46
47
    {}

    template <class Context, class QuadratureRule,
              class ElementMatrix, class RowNode, class ColNode, bool sameFE>
    void calculateElementMatrix(Context const& context,
                                QuadratureRule const& quad,
                                ElementMatrix& elementMatrix,
                                RowNode const& rowNode,
                                ColNode const& colNode,
                                std::integral_constant<bool, sameFE>,
                                std::false_type /*sameNode*/)
    {
48
49
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only.");
50
51
52
53
54
55
56
57
58
59
60
61
62
63

      auto const& rowLocalFE = rowNode.finiteElement();
      auto const& colLocalFE = colNode.finiteElement();

      test_exit( sameFE, "Not implemented: currently only the implementation for equal fespaces available" );

      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
        // Position of the current quadrature point in the reference element
        decltype(auto) local = context.position(quad[iq].position());

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

        // The multiplicative factor in the integral transformation formula
64
65
        double factor = context.integrationElement(quad[iq].position()) * quad[iq].weight();
        const auto exprValue = Super::coefficient(local);
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97

        // The gradients of the shape functions on the reference element
        std::vector<Dune::FieldMatrix<double,1,Context::dim> > referenceGradients;
        rowLocalFE.localBasis().evaluateJacobian(local, referenceGradients);

        // Compute the shape function gradients on the real element
        std::vector<Dune::FieldVector<double,Context::dow> > gradients(referenceGradients.size());

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

        for (std::size_t i = 0; i < rowLocalFE.size(); ++i) {
          const int local_i = rowNode.localIndex(i);
          for (std::size_t j = 0; j < colLocalFE.size(); ++j) {
            const int local_j = colNode.localIndex(j);
            elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
          }
        }
      }
    }


    template <class Context, class QuadratureRule,
              class ElementMatrix, class RowNode, class ColNode>
    void calculateElementMatrix(Context const& context,
                                QuadratureRule const& quad,
                                ElementMatrix& elementMatrix,
                                RowNode const& node,
                                ColNode const& /*colNode*/,
                                std::true_type /*sameFE*/,
                                std::true_type /*sameNode*/)
    {
98
99
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only.");
100

101
      using Category = ValueCategory_t<typename GridFct::Range>;
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
      calcElementMatrix(context, quad, elementMatrix, node, Category{});
    }

  protected:

    template <class Context, class QuadratureRule,
              class ElementMatrix, class RowNode>
    void calcElementMatrix(Context const& context,
                           QuadratureRule const& quad,
                           ElementMatrix& elementMatrix,
                           RowNode const& node,
                           tag::scalar)
    {
      auto const& localFE = node.finiteElement();
      std::vector<Dune::FieldMatrix<double,1,Context::dim> > referenceGradients;

      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
        // Position of the current quadrature point in the reference element
        decltype(auto) local = context.position(quad[iq].position());

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

        // The multiplicative factor in the integral transformation formula
126
        double factor = Super::coefficient(local) * context.integrationElement(quad[iq].position()) * quad[iq].weight();
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172

        // The gradients of the shape functions on the reference element
        localFE.localBasis().evaluateJacobian(local, referenceGradients);

        // Compute the shape function gradients on the real element
        std::vector<Dune::FieldVector<double,Context::dow> > gradients(referenceGradients.size());

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

        for (std::size_t i = 0; i < localFE.size(); ++i) {
          const int local_i = node.localIndex(i);

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

          for (std::size_t j = i+1; j < localFE.size(); ++j) {
            const int local_j = node.localIndex(j);
            double value = factor * (gradients[i] * gradients[j]);

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


    template <class Context, class QuadratureRule,
              class ElementMatrix, class RowNode>
    void calcElementMatrix(Context const& context,
                           QuadratureRule const& quad,
                           ElementMatrix& elementMatrix,
                           RowNode const& node,
                           tag::matrix)
    {
      auto const& localFE = node.finiteElement();
      std::vector<Dune::FieldMatrix<double,1,Context::dim> > referenceGradients;

      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
        // Position of the current quadrature point in the reference element
        decltype(auto) local = context.position(quad[iq].position());

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

        // The multiplicative factor in the integral transformation formula
173
174
        double factor = context.integrationElement(quad[iq].position()) * quad[iq].weight();
        const auto exprValue = Super::coefficient(local);
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213

        // The gradients of the shape functions on the reference element
        localFE.localBasis().evaluateJacobian(local, referenceGradients);

        // Compute the shape function gradients on the real element
        std::vector<Dune::FieldVector<double,Context::dow> > gradients(referenceGradients.size());

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

        for (std::size_t i = 0; i < localFE.size(); ++i) {
          const int local_i = node.localIndex(i);
          for (std::size_t j = 0; j < localFE.size(); ++j) {
            const int local_j = node.localIndex(j);
            elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
          }
        }
      }
    }

  protected:

    template <int dow>
    double eval(Dune::FieldVector<double,1> const& scalar, double factor,
                Dune::FieldVector<double,dow> const& grad_test,
                Dune::FieldVector<double,dow> const& grad_trial) const
    {
      return (scalar * factor) * (grad_test * grad_trial);
    }

    template <int dow>
    double eval(Dune::FieldMatrix<double,dow,dow> const& M, double factor,
                Dune::FieldVector<double,dow> const& grad_test,
                Dune::FieldVector<double,dow> const& grad_trial) const
    {
      return factor * (grad_test * (M * grad_trial));
    }
  };

214
215
  /** @} **/

216
} // end namespace AMDiS