ConvectionDiffusionOperator.hpp 10.2 KB
Newer Older
1
2
3
4
#pragma once

#include <type_traits>

5
#include <amdis/utility/LocalBasisCache.hpp>
6
7
#include <amdis/LocalOperator.hpp>
#include <amdis/Output.hpp>
8
9
#include <amdis/common/TypeTraits.hpp>
#include <amdis/common/StaticSize.hpp>
10
11
12
13
14
15
16
17

namespace AMDiS
{
  /**
   * \addtogroup operators
   * @{
   **/

18
  /// convection-diffusion operator, see \ref convectionDiffusion
19
20
  /// <A*grad(u),grad(v)> - <b*u, grad(v)> + <c*u, v> = <f, v> (conserving) or
  /// <A*grad(u),grad(v)> + <b*grad(u), v> + <c*u, v> = <f, v> (non conserving)
21
  template <class LocalContext, class GridFctA, class GridFctB, class GridFctC, class GridFctF, bool conserving = true>
22
  class ConvectionDiffusionOperator
23
24
      : public LocalOperator<ConvectionDiffusionOperator<LocalContext, GridFctA, GridFctB, GridFctC, GridFctF, conserving>,
                             LocalContext>
25
  {
26
27
28
29
    static const int dow = ContextGeometry<LocalContext>::dow;

    using A_range = typename GridFctA::Range;
    static_assert( Size_v<A_range> == 1 || (Rows_v<A_range> == dow && Cols_v<A_range> == dow),
30
      "Expression A(x) must be of scalar or matrix type." );
31
32
    using b_range = typename GridFctB::Range;
    static_assert( Size_v<b_range> == 1 || Size_v<b_range> == dow,
33
      "Expression b(x) must be of scalar or vector type." );
34
35
    using c_range = typename GridFctC::Range;
    static_assert( Size_v<c_range> == 1,
36
      "Expression c(x) must be of scalar type." );
37
38
    using f_range = typename GridFctF::Range;
    static_assert( Size_v<f_range> == 1,
39
40
41
42
      "Expression f(x) must be of scalar type." );

  public:
    ConvectionDiffusionOperator(GridFctA const& gridFctA, GridFctB const& gridFctB,
43
                                GridFctC const& gridFctC, GridFctF const& gridFctF)
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
      : gridFctA_(gridFctA)
      , gridFctB_(gridFctB)
      , gridFctC_(gridFctC)
      , gridFctF_(gridFctF)
    {}

    template <class Context, class RowNode, class ColNode, class ElementMatrix>
    void getElementMatrix(Context const& context,
                          RowNode const& rowNode, ColNode const& colNode,
                          ElementMatrix& elementMatrix)
    {
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only." );

      static_assert(std::is_same<FiniteElementType_t<RowNode>, FiniteElementType_t<ColNode>>{},
        "Galerkin operator requires equal finite elements for test and trial space." );

61
      using RangeFieldType = typename NodeQuadCache<RowNode>::Traits::RangeFieldType;
62
63
64
65
66
67

      auto localFctA = localFunction(gridFctA_); localFctA.bind(context.element());
      auto localFctB = localFunction(gridFctB_); localFctB.bind(context.element());
      auto localFctC = localFunction(gridFctC_); localFctC.bind(context.element());

      auto const& localFE = rowNode.finiteElement();
68
      std::size_t numLocalFe = localFE.size();
69
70
71
72
73

      int quadDeg = std::max({this->getDegree(2,coeffOrder(localFctA),rowNode,colNode),
                              this->getDegree(1,coeffOrder(localFctB),rowNode,colNode),
                              this->getDegree(0,coeffOrder(localFctC),rowNode,colNode)});

74
      using QuadratureRules = Dune::QuadratureRules<typename Context::Geometry::ctype, Context::LocalContext::mydimension>;
75
76
      auto const& quad = QuadratureRules::rule(context.type(), quadDeg);

77
78
79
      NodeQuadCache<RowNode> cache(localFE.localBasis());
      auto const& shapeGradientsCache = cache.evaluateJacobianAtQP(context, quad);
      auto const& shapeValuesCache = cache.evaluateFunctionAtQP(context, quad);
80
81

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

        // 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
89
        const auto factor = context.integrationElement(quad[iq].position()) * quad[iq].weight();
90
91

        // the values of the shape functions on the reference element at the quadrature point
92
        auto const& shapeValues = shapeValuesCache[iq];
93
94

        // The gradients of the shape functions on the reference element
95
        auto const& shapeGradients = shapeGradientsCache[iq];
96
97

        // Compute the shape function gradients on the real element
98
        using WorldVector = FieldVector<RangeFieldType,Context::dow>;
99
100
        thread_local std::vector<WorldVector> gradients;
        gradients.resize(shapeGradients.size());
101
102
103
104
105

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

        const auto A = localFctA(local);
106
        WorldVector b = makeB(localFctB(local));
107
108
109
        const auto c = localFctC(local);

        IF_CONSTEXPR(conserving) {
110
          WorldVector gradAi;
111
          for (std::size_t i = 0; i < numLocalFe; ++i) {
112
            const auto local_i = rowNode.localIndex(i);
113
            gradAi = A * gradients[i];
114
            auto gradBi = b * gradients[i];
115
            for (std::size_t j = 0; j < numLocalFe; ++j) {
116
              const auto local_j = colNode.localIndex(j);
117
118
119
120
              elementMatrix[local_i][local_j] += (dot(gradAi, gradients[j]) + (c * shapeValues[i] - gradBi) * shapeValues[j]) * factor;
            }
          }
        } else {
121
          WorldVector grad_i;
122
          for (std::size_t i = 0; i < numLocalFe; ++i) {
123
            const auto local_i = rowNode.localIndex(i);
124
125
            grad_i = A * gradients[i];
            grad_i+= b * shapeValues[i];
126
            for (std::size_t j = 0; j < numLocalFe; ++j) {
127
              const auto local_j = colNode.localIndex(j);
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
              elementMatrix[local_i][local_j] += (dot(grad_i, gradients[j]) + c * shapeValues[i] * shapeValues[j]) * factor;
            }
          }
        }
      }

      localFctA.unbind();
      localFctB.unbind();
      localFctC.unbind();
    }


    template <class Context, class Node, class ElementVector>
    void getElementVector(Context const& context,
                          Node const& node,
                          ElementVector& elementVector)
    {
145
      static_assert(Node::isLeaf,
146
147
148
149
150
        "Operator can be applied to Leaf-Nodes only." );

      auto localFctF = localFunction(gridFctF_); localFctF.bind(context.element());

      auto const& localFE = node.finiteElement();
151
      std::size_t numLocalFe = localFE.size();
152
153
154

      int quad_order = this->getDegree(0,coeffOrder(localFctF),node);

155
      using QuadratureRules = Dune::QuadratureRules<typename Context::Geometry::ctype, Context::LocalContext::dimension>;
156
157
      auto const& quad = QuadratureRules::rule(context.type(), quad_order);

158
159
      NodeQuadCache<Node> cache(localFE.localBasis());
      auto const& shapeValuesCache = cache.evaluateFunctionAtQP(context, quad);
160

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

        // the values of the shape functions on the reference element at the quadrature point
166
        auto const& shapeValues = shapeValuesCache[iq];
167

168
169
        // The multiplicative factor in the integral transformation formula
        const auto factor = localFctF(local) * context.integrationElement(quad[iq].position()) * quad[iq].weight();
170

171
        for (std::size_t i = 0; i < numLocalFe; ++i) {
172
          const auto local_i = node.localIndex(i);
173
          elementVector[local_i] += shapeValues[i] * factor;
174
175
176
177
178
179
180
181
182
183
184
185
186
187
        }
      }

      localFctF.unbind();
    }

  private:

    template <class LF>
    using HasLocalFunctionOrder = decltype( order(std::declval<LF>()) );

    template <class LocalFct>
    int coeffOrder(LocalFct const& localFct)
    {
188
      using Concept = Dune::Std::is_detected<HasLocalFunctionOrder, LocalFct>;
189
190
191
192
193
      return Dune::Hybrid::ifElse(Concept{},
        [&](auto id) { return order(id(localFct)); },
        [] (auto)    { return 0; });
    }

194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
    template <class T, int N>
    static FieldVector<T,dow> makeB(FieldVector<T,N> const& b) { return b; }

    template <class T, int N>
    static FieldVector<T,dow> makeB(FieldVector<T,N>&& b) { return std::move(b); }

    template <class T>
    static FieldVector<T,dow> makeB(FieldVector<T,1> const& b) { return {T(b)}; }

    template <class T>
    static FieldVector<T,dow> makeB(FieldVector<T,1>&& b) { return {T(b)}; }

    template <class T, std::enable_if_t<std::is_arithmetic<T>::value, int> = 0>
    static FieldVector<T,dow> makeB(T b) { return {b}; }

209
  private:
210
211
212
213
    GridFctA gridFctA_;
    GridFctB gridFctB_;
    GridFctC gridFctC_;
    GridFctF gridFctF_;
214
215
  };

216
  template <class PreGridFctA, class PreGridFctB, class PreGridFctC, class PreGridFctF, class c>
217
218
  struct PreConvectionDiffusionOperator
  {
219
    static constexpr bool conserving = c::value;
220
221
222
223
224
    PreGridFctA gridFctA;
    PreGridFctB gridFctB;
    PreGridFctC gridFctC;
    PreGridFctF gridFctF;
  };
225

226
227
228
  template <class PreGridFctA, class PreGridFctB, class PreGridFctC, class PreGridFctF, bool conserving = true>
  auto convectionDiffusion(PreGridFctA const& gridFctA, PreGridFctB const& gridFctB,
                           PreGridFctC const& gridFctC, PreGridFctF const& gridFctF,
229
230
                           bool_t<conserving> = {})
  {
231
    using Pre = PreConvectionDiffusionOperator<PreGridFctA, PreGridFctB, PreGridFctC, PreGridFctF, bool_t<conserving>>;
232
233
234
    return Pre{gridFctA, gridFctB, gridFctC, gridFctF};
  }

235
236
  template <class Context, class... T, class GridView>
  auto makeLocalOperator(PreConvectionDiffusionOperator<T...> pre, GridView const& gridView)
237
  {
238
    using Pre = PreConvectionDiffusionOperator<T...>;
239
240
241
242
    auto gridFctA = makeGridFunction(std::move(pre.gridFctA), gridView);
    auto gridFctB = makeGridFunction(std::move(pre.gridFctB), gridView);
    auto gridFctC = makeGridFunction(std::move(pre.gridFctC), gridView);
    auto gridFctF = makeGridFunction(std::move(pre.gridFctF), gridView);
243

244
    using GridFctOp = ConvectionDiffusionOperator<Context,
245
      decltype(gridFctA), decltype(gridFctB), decltype(gridFctC), decltype(gridFctF), Pre::conserving>;
246

247
    GridFctOp localOperator{std::move(gridFctA), std::move(gridFctB), std::move(gridFctC), std::move(gridFctF)};
248
    return localOperator;
249
250
251
252
253
  }

  /** @} **/

} // end namespace AMDiS