Liebe Gitlab-Nutzerin, lieber Gitlab-Nutzer,
es ist nun möglich sich mittels des ZIH-Logins/LDAP an unserem Dienst anzumelden. Die Konten der externen Nutzer:innen sind über den Reiter "Standard" erreichbar.
Die Administratoren


Dear Gitlab user,
it is now possible to log in to our service using the ZIH login/LDAP. The accounts of external users can be accessed via the "Standard" tab.
The administrators

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

#include <type_traits>

5
#include <amdis/LocalBasisCache.hpp>
6 7 8
#include <amdis/LocalOperator.hpp>
#include <amdis/Output.hpp>
#include <amdis/common/Utility.hpp>
9
#include <amdis/common/Size.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 61 62 63 64 65 66 67 68
      : 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." );

      using LocalBasisType = typename FiniteElementType_t<RowNode>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;

      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();
69
      std::size_t numLocalFe = localFE.size();
70 71 72 73 74

      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)});

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

78 79 80 81 82
      LocalBasisCache<LocalBasisType> cache(localFE.localBasis());
      auto const& shapeGradientsCache = cache.evaluateJacobianAtQp(context, quad);
      auto const& shapeValuesCache = cache.evaluateFunctionAtQp(context, quad);

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

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

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

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

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

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

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

        IF_CONSTEXPR(conserving) {
111
          WorldVector gradAi, gradBi;
112
          for (std::size_t i = 0; i < numLocalFe; ++i) {
113
            const int local_i = rowNode.localIndex(i);
114 115
            gradAi = A * gradients[i];
            gradBi = b * gradients[i];
116
            for (std::size_t j = 0; j < numLocalFe; ++j) {
117 118 119 120 121
              const int local_j = colNode.localIndex(j);
              elementMatrix[local_i][local_j] += (dot(gradAi, gradients[j]) + (c * shapeValues[i] - gradBi) * shapeValues[j]) * factor;
            }
          }
        } else {
122
          WorldVector grad_i;
123
          for (std::size_t i = 0; i < numLocalFe; ++i) {
124
            const int local_i = rowNode.localIndex(i);
125 126
            grad_i = A * gradients[i];
            grad_i+= b * shapeValues[i];
127
            for (std::size_t j = 0; j < numLocalFe; ++j) {
128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145
              const int local_j = colNode.localIndex(j);
              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)
    {
146
      static_assert(Node::isLeaf,
147 148
        "Operator can be applied to Leaf-Nodes only." );

149
      using LocalBasisType = typename FiniteElementType_t<Node>::Traits::LocalBasisType;
150 151 152 153

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

      auto const& localFE = node.finiteElement();
154
      std::size_t numLocalFe = localFE.size();
155 156 157

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

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

161 162
      LocalBasisCache<LocalBasisType> cache(localFE.localBasis());
      auto const& shapeValuesCache = cache.evaluateFunctionAtQp(context, quad);
163

164 165 166
      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());
167 168

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

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

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

      localFctF.unbind();
    }

  private:

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

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

  private:
198 199 200 201
    GridFctA gridFctA_;
    GridFctB gridFctB_;
    GridFctC gridFctC_;
    GridFctF gridFctF_;
202 203
  };

204
  template <class PreGridFctA, class PreGridFctB, class PreGridFctC, class PreGridFctF, class c>
205 206
  struct PreConvectionDiffusionOperator
  {
207
    static constexpr bool conserving = c::value;
208 209 210 211 212
    PreGridFctA gridFctA;
    PreGridFctB gridFctB;
    PreGridFctC gridFctC;
    PreGridFctF gridFctF;
  };
213

214 215 216
  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,
217 218
                           bool_t<conserving> = {})
  {
219
    using Pre = PreConvectionDiffusionOperator<PreGridFctA, PreGridFctB, PreGridFctC, PreGridFctF, bool_t<conserving>>;
220 221 222
    return Pre{gridFctA, gridFctB, gridFctC, gridFctF};
  }

223 224
  template <class Context, class... T, class GridView>
  auto makeLocalOperator(PreConvectionDiffusionOperator<T...> pre, GridView const& gridView)
225
  {
226
    using Pre = PreConvectionDiffusionOperator<T...>;
227 228 229 230
    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);
231

232
    using GridFctOp = ConvectionDiffusionOperator<Context,
233
      decltype(gridFctA), decltype(gridFctB), decltype(gridFctC), decltype(gridFctF), Pre::conserving>;
234

235
    GridFctOp localOperator{std::move(gridFctA), std::move(gridFctB), std::move(gridFctC), std::move(gridFctF)};
236
    return localOperator;
237 238 239 240 241
  }

  /** @} **/

} // end namespace AMDiS