GridFunctionOperator.hpp 9.16 KB
Newer Older
1
2
3
#pragma once

#include <cassert>
4
#include <type_traits>
5

6
#include <amdis/GridFunctions.hpp>
7
#include <amdis/LocalOperator.hpp>
8
#include <amdis/common/Transposed.hpp>
9
10
11
#include <amdis/common/TypeTraits.hpp>
#include <amdis/typetree/FiniteElementType.hpp>
#include <amdis/utility/QuadratureFactory.hpp>
12
13
14

namespace AMDiS
{
15
16
17
18
19
  /**
   * \addtogroup operators
   * @{
   **/

20
  /// \brief The main implementation of an CRTP-base class for operators using a grid-function
21
  /// coefficient to be used in an \ref Assembler.
22
23
24
25
26
27
28
  /**
   * An Operator that takes a GridFunction as coefficient.
   * Provides quadrature rules and handles the evaluation of the GridFunction at
   * local coordinates.
   *
   * The class is specialized, by deriving from it, in \ref GridFunctionOperator.
   *
29
   * \tparam Derived      The class derived from GridFunctionOperatorBase
30
   * \tparam LocalContext The Element or Intersection type
31
32
33
34
   * \tparam GridFunction The GridFunction, a LocalFunction is created from, and
   *                      that is evaluated at quadrature points.
   *
   * **Requirements:**
35
   * - `LocalContext` models either Entity (of codim 0) or Intersection
36
37
   * - `GridFunction` models the \ref Concepts::GridFunction
   **/
38
  template <class Derived, class LocalContext, class GridFunction>
39
  class GridFunctionOperatorBase
40
      : public LocalOperator<Derived, LocalContext>
41
  {
42
43
44
45
46
47
48
    template <class D, class LC>
    friend class LocalOperator;

    using ContextType = Impl::ContextTypes<LocalContext>;
    using Super = LocalOperator<Derived, LocalContext>;

  private:
49
    /// The type of the localFunction associated with the GridFunction
50
    using LocalFunction = decltype(localFunction(std::declval<GridFunction>()));
51

52
    /// The Codim=0 entity of the grid, the localFunction can be bound to
53
54
    using Element = typename ContextType::Entity;

55
    /// The geometry-type of the grid element
56
57
    using Geometry = typename Element::Geometry;

58
    /// The type of the local coordinates in the \ref Element
59
60
    using LocalCoordinate = typename GridFunction::EntitySet::LocalCoordinate;

61
    /// A factory for QuadratureRules that incooperate the order of the LocalFunction
62
    using QuadFactory = QuadratureFactory<typename Geometry::ctype, LocalContext::mydimension, LocalFunction>;
63

64
  public:
65
    /// \brief Constructor. Stores a copy of `gridFct`.
66
    /**
67
68
69
     * A GridFunctionOperator takes a gridFunction and the
     * differentiation order of the operator, to calculate the
     * quadrature degree in \ref getDegree.
70
     **/
71
72
    template <class GF>
    GridFunctionOperatorBase(GF&& gridFct, int termOrder)
73
      : gridFct_(FWD(gridFct))
74
      , termOrder_(termOrder)
75
76
    {}

77
    /// Create a quadrature factory from a PreQuadratureFactory, e.g. class derived from \ref QuadratureFactory
78
    template <class PreQuadFactory>
79
    void setQuadFactory(PreQuadFactory&& pre)
80
    {
81
      using ctype = typename Geometry::ctype;
82
      quadFactory_.emplace(
83
        makeQuadratureFactory<ctype, LocalContext::mydimension, LocalFunction>(FWD(pre)));
84
85
86
    }

  protected:
87
    /// Return expression value at LocalCoordinates
88
    auto coefficient(LocalCoordinate const& local) const
89
    {
90
      assert( this->bound_ );
91
      return (*localFct_)(local);
92
93
    }

94
    /// Create a quadrature rule using the \ref QuadratureFactory by calculating the
95
96
    /// quadrature order necessary to integrate the (bi)linear-form accurately.
    template <class... Nodes>
97
    auto const& getQuadratureRule(Dune::GeometryType type, Nodes const&... nodes) const
98
    {
99
100
101
      assert( bool(quadFactory_) );
      int quadDegree = this->getDegree(termOrder_, quadFactory_->order(), nodes...);
      return quadFactory_->rule(type, quadDegree);
102
103
    }

104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
  private:
    /// \brief Binds operator to `element` and `geometry`.
    /**
     * Binding an operator to the currently visited element in grid traversal.
     * Since all operators need geometry information, the `Element::Geometry`
     * object `geometry` is created once, during grid traversal, and passed in.
     *
     * By default, it binds the \ref localFct_ to the `element` and the Quadrature
     * factory to the localFunction.
     **/
    void bind_impl(Element const& element, Geometry const& geometry)
    {
      assert( bool(quadFactory_) );
      localFct_.emplace(localFunction(gridFct_));
      localFct_->bind(element);
      quadFactory_->bind(localFct_.value());
    }

    /// Unbinds operator from element.
    void unbind_impl()
    {
      localFct_->unbind();
      localFct_.reset();
    }

129
  private:
130
    /// The gridFunction to be used within the operator
131
    GridFunction gridFct_;
132
133

    /// localFunction associated with gridFunction. Mus be updated whenever gridFunction changes.
134
    Dune::Std::optional<LocalFunction> localFct_;
135

136
    /// Assign each element type a quadrature rule
137
    Dune::Std::optional<QuadFactory> quadFactory_;
138

139
140
    /// the derivative order of this operator (in {0, 1, 2})
    int termOrder_ = 0;
141
142
143
  };


144
145
  /// \brief The transposed operator, implemented in term of its transposed by
  /// calling \ref getElementMatrix with inverted arguments.
146
147
  template <class Derived, class Transposed>
  class GridFunctionOperatorTransposed
148
      : public LocalOperator<Derived, typename Transposed::LocalContext>
149
  {
150
151
152
    template <class D, class LC>
    friend class LocalOperator;

153
154
155
    template <class T, class... Args>
    using Constructable = decltype( new T(std::declval<Args>()...) );

156
  public:
157
158
159
    template <class... Args,
      std::enable_if_t<Dune::Std::is_detected<Constructable, Transposed, Args...>::value, int> = 0>
    GridFunctionOperatorTransposed(Args&&... args)
160
      : transposedOp_(FWD(args)...)
161
162
    {}

163
164
165
166
    /// Redirects the setQuadFactory call top the transposed operator
    template <class PreQuadFactory>
    void setQuadFactory(PreQuadFactory&& pre)
    {
167
      transposedOp_.setQuadFactory(FWD(pre));
168
169
170
    }

  private:
171
    /// Redirects the bind call top the transposed operator
172
173
    template <class Element, class Geometry>
    void bind_impl(Element const& element, Geometry const& geometry)
174
    {
175
      transposedOp_.bind(element, geometry);
176
177
    }

178
    /// Redirects the unbind call top the transposed operator
179
    void unbind_impl()
180
    {
181
      transposedOp_.unbind();
182
183
    }

184
    /// Apply the assembling to the transposed elementMatrix with interchanged row-/colNode
185
186
187
188
189
190
    template <class Context, class RowNode, class ColNode, class ElementMatrix>
    void getElementMatrix(Context const& context,
                          RowNode const& rowNode,
                          ColNode const& colNode,
                          ElementMatrix& elementMatrix)
    {
191
      auto elementMatrixTransposed = transposed(elementMatrix);
192
193
194
      transposedOp_.getElementMatrix(
        context, colNode, rowNode, elementMatrixTransposed);
    }
195

196
197
198
  private:
    Transposed transposedOp_;
  };
199

200

201
#ifndef DOXYGEN
202
203
  template <class Tag, class PreGridFct, class PreQuadFactory>
  struct PreGridFunctionOperator
204
205
  {
    Tag tag;
206
207
    PreGridFct expr;
    PreQuadFactory quadFactory;
208
  };
209
#endif
210

211
212
213
  /// Store tag and expression into a \ref PreGridFunctionOperator to create a \ref GridFunctionOperator
  template <class Tag, class Expr, class... QuadratureArgs>
  auto makeOperator(Tag tag, Expr&& expr, QuadratureArgs&&... args)
214
  {
215
216
    auto pqf = makePreQuadratureFactory(FWD(args)...);
    return PreGridFunctionOperator<Tag, TYPEOF(expr), TYPEOF(pqf)>{tag, FWD(expr), std::move(pqf)};
217
218
  }

219
220
221
222
  /** @} **/

#ifndef DOXYGEN

223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
  /// The base-template for GridFunctionOperators
  /**
   * An operator can specialize this class, by deriving from \ref GridFunctionOperatorBase.
   * With the generic function \ref makeLocalOperator, an instance is created. To
   * distinguisch different GridFunction operators, a tag can be provided that has no
   * other effect.
   *
   * \tparam Tag  An Identifier for this GridFunctionOperator
   * \tparam LocalContext  An Element or Intersection the operator is evaluated on
   * \tparam GridFct  A GridFunction evaluated in local coordinates on the bound element
   **/
  template <class Tag, class LocalContext, class GridFct>
  class GridFunctionOperator
      : public GridFunctionOperatorBase<GridFunctionOperator<Tag, LocalContext, GridFct>, LocalContext, GridFct>
  {};

239
240
  template <class Context, class Tag, class GF, class QF>
  auto makeGridFunctionOperator(Tag tag, GF&& gf, QF&& qf)
241
  {
242
243
    GridFunctionOperator<Tag, Context, TYPEOF(gf)> gfo{tag, FWD(gf)};
    gfo.setQuadFactory(FWD(qf));
244
    return gfo;
245
246
  }

247

248
  /// Generate an \ref GridFunctionOperator from a PreOperator (tag, expr).
249
  /// @{
250
251
  template <class Context, class... Args, class GridView>
  auto makeLocalOperator(PreGridFunctionOperator<Args...> op, GridView const& gridView)
252
  {
253
254
    auto gf = makeGridFunction(std::move(op.expr), gridView);
    return makeGridFunctionOperator<Context>(op.tag, std::move(gf), std::move(op.quadFactory));
255
256
  }

257
258
  template <class Context, class... Args, class GridView>
  auto makeLocalOperator(std::reference_wrapper<PreGridFunctionOperator<Args...>> op_ref, GridView const& gridView)
259
  {
260
261
262
    PreGridFunctionOperator<Args...> const& op = op_ref;
    auto gf = makeGridFunction(std::ref(op.expr), gridView);
    return makeGridFunctionOperator<Context>(op.tag, std::move(gf), op.quadFactory);
263
  }
264
  /// @}
265

266
267
#endif // DOXYGEN

268
} // end namespace AMDiS