AnalyticGridFunction.hpp 5.81 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
3
4
5
#pragma once

#include <type_traits>

#include <dune/common/std/optional.hh>
6
#include <dune/functions/common/signature.hh>
Praetorius, Simon's avatar
Praetorius, Simon committed
7
8
#include <dune/functions/gridfunctions/gridviewentityset.hh>

9
10
#include <amdis/Operations.hpp>
#include <amdis/gridfunctions/GridFunctionConcepts.hpp>
Praetorius, Simon's avatar
Praetorius, Simon committed
11
12
13

namespace AMDiS
{
14
#ifndef DOXYGEN
15
16
  template <class Signature, class LocalContext, class Function>
  class AnalyticLocalFunction;
17
#endif
Praetorius, Simon's avatar
Praetorius, Simon committed
18

19
20
  template <class R, class D, class LocalContext, class Function>
  class AnalyticLocalFunction<R(D), LocalContext, Function>
Praetorius, Simon's avatar
Praetorius, Simon committed
21
22
  {
  public:
23
24
    using Range = R;
    using Domain = D; // LocalDomain
Praetorius, Simon's avatar
Praetorius, Simon committed
25

26
    using Geometry = typename LocalContext::Geometry;
Praetorius, Simon's avatar
Praetorius, Simon committed
27

28
29
    enum { hasDerivative = true };

Praetorius, Simon's avatar
Praetorius, Simon committed
30
  public:
31
    AnalyticLocalFunction(Function const& fct)
32
      : fct_{fct}
33
34
35
36
37
38
    {}

    void bind(LocalContext const& element)
    {
      geometry_.emplace(element.geometry());
    }
Praetorius, Simon's avatar
Praetorius, Simon committed
39

40
    void unbind()
Praetorius, Simon's avatar
Praetorius, Simon committed
41
    {
42
43
      geometry_.reset();
    }
Praetorius, Simon's avatar
Praetorius, Simon committed
44

45
46
    Range operator()(Domain const& local) const
    {
47
      assert( bool(geometry_) );
48
49
50
51
52
53
54
55
56
57
58
59
60
      return fct_(geometry_.value().global(local));
    }

    Function const& fct() const
    {
      return fct_;
    }

  private:
    Function fct_;
    Dune::Std::optional<Geometry> geometry_;
  };

61
62
  /// \brief Return the polynomial order of the function f, if a free function
  /// order(f) exists, otherwise return a default value. \relates AnalyticLocalFunction
63
64
  template <class Sig, class LocalContext, class F,
    REQUIRES(Concepts::HasOrder<F>)>
65
66
  int order(AnalyticLocalFunction<Sig,LocalContext,F> const& lf)
  {
67
    return order(lf.fct(),1);
68
  }
69
70


71
72
73
74
75
76
77
78
79
80
81
82
83
84
  template <class R, class D, class LocalContext, class F>
  auto derivative(AnalyticLocalFunction<R(D),LocalContext,F> const& lf)
  {
    static_assert(Concepts::HasPartial<F>,
      "No partial(_0,...) defined for Functor F of AnalyticLocalFunction.");

    auto df = partial(lf.fct(), index_<0>);

    using RawSignature = typename Dune::Functions::SignatureTraits<R(D)>::RawSignature;
    using DerivativeSignature = typename Dune::Functions::DefaultDerivativeTraits<RawSignature>::Range(D);
    return AnalyticLocalFunction<DerivativeSignature,LocalContext,decltype(df)>{df};
  }


85
86
87
88
89
90
91
92
93
94
95
96
  /// \class AnalyticGridFunction
  /// \brief A Gridfunction that evaluates a function with global coordinates.
  /**
   * \ingroup GridFunctions
   * Implements a GridFunction that wraps a functor around global coordinates.
   *
   * \tparam Function The callable `f=f(x)` with `x in Domain == GlobalCoordinates`
   * \tparam GridView A GridView that models `Dune::GridViewConcept`.
   *
   * **Requirements:**
   * - Function models \ref Concepts::Callable<Domain>
   **/
97
98
99
100
101
102
103
  template <class Function, class GridView>
  class AnalyticGridFunction
  {
  public:
    using EntitySet = Dune::Functions::GridViewEntitySet<GridView, 0>;
    using Domain = typename EntitySet::GlobalCoordinate;
    using Range = std::decay_t<std::result_of_t<Function(Domain)>>;
Praetorius, Simon's avatar
Praetorius, Simon committed
104

105
106
    enum { hasDerivative = true };

107
108
109
  private:
    using Element = typename EntitySet::Element;
    using LocalDomain = typename EntitySet::LocalCoordinate;
110
    using LocalFunction = AnalyticLocalFunction<Range(LocalDomain), Element, Function>;
Praetorius, Simon's avatar
Praetorius, Simon committed
111

112
  public:
113
    /// \brief Constructor. Stores the function `fct` and creates an `EntitySet`.
Praetorius, Simon's avatar
Praetorius, Simon committed
114
    AnalyticGridFunction(Function const& fct, GridView const& gridView)
115
116
      : fct_{fct}
      , entitySet_{gridView}
Praetorius, Simon's avatar
Praetorius, Simon committed
117
118
    {}

119
    /// \brief Return the evaluated functor at global coordinates
Praetorius, Simon's avatar
Praetorius, Simon committed
120
121
122
123
124
    Range operator()(Domain const& x) const
    {
      return fct_(x);
    }

125
    /// \brief Return the LocalFunction of the AnalyticGridFunction.
126
    LocalFunction localFunction() const
Praetorius, Simon's avatar
Praetorius, Simon committed
127
    {
128
      return {fct_};
Praetorius, Simon's avatar
Praetorius, Simon committed
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
    }

    EntitySet const& entitySet() const
    {
      return entitySet_;
    }

    Function const& fct() const { return fct_; }

  private:
    Function fct_;

    EntitySet entitySet_;
  };


145
146
147
148
149
150
  template <class F, class GV>
  auto localFunction(AnalyticGridFunction<F,GV> const& gf)
  {
    return gf.localFunction();
  }

151
152
153
154
155
156
157
158
159
160
161
162
163
  /// \brief Return a GridFunction representing the derivative of a functor.
  // [expects: Functor f has free function derivative(f)]
  template <class F, class GV>
  auto derivative(AnalyticGridFunction<F,GV> const& gf)
  {
    static_assert(Concepts::HasPartial<F>,
      "No partial(_0,...) defined for Functor of AnalyticLocalFunction.");

    auto df = partial(gf.fct(), index_<0>);
    return AnalyticGridFunction<decltype(df), GV>{df, gf.entitySet().gridView()};
  }


164
165
#ifndef DOXYGEN
  // A pre-GridFunction that just stores the function
Praetorius, Simon's avatar
Praetorius, Simon committed
166
167
168
  template <class Function>
  struct AnalyticPreGridFunction
  {
169
170
171
172
173
174
175
176
177
178
179
    using Self = AnalyticPreGridFunction;

    struct Creator
    {
      template <class GridView>
      static auto create(Self const& self, GridView const& gridView)
      {
        return AnalyticGridFunction<Function, GridView>{self.fct_, gridView};
      }
    };

Praetorius, Simon's avatar
Praetorius, Simon committed
180
181
182
    Function fct_;
  };

183
  namespace Traits
Praetorius, Simon's avatar
Praetorius, Simon committed
184
185
186
187
  {
    template <class Functor>
    struct IsPreGridFunction<AnalyticPreGridFunction<Functor>>
      : std::true_type {};
188
  }
189
#endif
190

Praetorius, Simon's avatar
Praetorius, Simon committed
191

192
193
194
195
196
197
198
199
200
  /// \fn evalAtQP \brief Generator function for AnalyticGridFunction. \relates AnalyticGridfunction
  /**
   * \ingroup GridFunctions
   * Evaluate a functor at global coordinates. See \ref AnalyticGridFunction.
   *
   * **Examples:**
   * - `evalAtQP([](Dune::FieldVector<double, 2> const& x) { return x[0]; })`
   * - `evalAtQP(Operation::TwoNorm{})`
   **/
Praetorius, Simon's avatar
Praetorius, Simon committed
201
202
  template <class Function,
    REQUIRES(Concepts::CallableDomain<Function>)>
203
  auto evalAtQP(Function const& f)
Praetorius, Simon's avatar
Praetorius, Simon committed
204
205
206
207
208
  {
    return AnalyticPreGridFunction<Function>{f};
  }


209
210
211
212
213
214
215
216
217
218
  template <class Function>
  struct GridFunctionCreator<Function, std::enable_if_t<Concepts::CallableDomain<Function>>>
  {
    template <class GridView>
    static auto create(Function const& fct, GridView const& gridView)
    {
      return AnalyticGridFunction<Function, GridView>{fct, gridView};
    }
  };

Praetorius, Simon's avatar
Praetorius, Simon committed
219
} // end namespace AMDiS