DiscreteLocalFunction.inc.hpp 11.9 KB
Newer Older
1
2
#pragma once

Praetorius, Simon's avatar
Praetorius, Simon committed
3
#include <amdis/Output.hpp>
4
5
#include <amdis/common/DerivativeTraits.hpp>
#include <amdis/common/FieldMatVec.hpp>
Praetorius, Simon's avatar
Praetorius, Simon committed
6
#include <amdis/functions/NodeIndices.hpp>
7
8
9
10
11
12
#include <amdis/utility/LocalBasisCache.hpp>
#include <amdis/utility/LocalToGlobalAdapter.hpp>

#include <dune/common/ftraits.hh>

namespace AMDiS {
Praetorius, Simon's avatar
Praetorius, Simon committed
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
namespace Impl {

// specialization of Coeff has gather method
template <class Coeff, class LocalView, class LocalCoeff,
  class = decltype(std::declval<Coeff>().gather(std::declval<LocalView>(), std::declval<LocalCoeff&>()))>
void gather(Coeff const& coeff, LocalView const& localView, LocalCoeff& localCoeff, Dune::PriorityTag<2>)
{
  coeff.gather(localView, localCoeff);
}

// fallback implementation
template <class Coeff, class LocalView, class LocalCoeff>
void gather(Coeff const& coeff, LocalView const& localView, LocalCoeff& localCoeff, Dune::PriorityTag<1>)
{
  localCoeff.resize(localView.size());
  auto it = localCoeff.begin();
  for (auto const& idx : nodeIndices(localView))
    *it++ = coeff[idx];
}

} // end namespace Impl

35

Praetorius, Simon's avatar
Praetorius, Simon committed
36
37
template <class Coeff, class GB, class TP>
class DiscreteFunction<Coeff const,GB,TP>::LocalFunction
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
{
public:
  using Domain = typename EntitySet::LocalCoordinate;
  using Range = typename DiscreteFunction::Range;

  enum { hasDerivative = true };

private:
  using LocalView = typename GlobalBasis::LocalView;
  using Element = typename EntitySet::Element;
  using Geometry = typename Element::Geometry;

public:
  /// Constructor. Stores a copy of the DiscreteFunction.
  LocalFunction(DiscreteFunction const& globalFunction)
    : globalFunction_(globalFunction)
Praetorius, Simon's avatar
Praetorius, Simon committed
54
    , localView_(globalFunction_.basis().localView())
55
    , subTree_(&Dune::TypeTree::child(localView_.tree(), globalFunction_.treePath()))
56
57
58
59
60
  {}

  /// Copy constructor.
  LocalFunction(LocalFunction const& other)
    : globalFunction_(other.globalFunction_)
Praetorius, Simon's avatar
Praetorius, Simon committed
61
    , localView_(globalFunction_.basis().localView())
62
    , subTree_(&Dune::TypeTree::child(localView_.tree(), globalFunction_.treePath()))
63
64
65
66
67
68
  {}

  /// \brief Bind the LocalView to the element
  void bind(Element const& element)
  {
    localView_.bind(element);
Praetorius, Simon's avatar
Praetorius, Simon committed
69
    Impl::gather(globalFunction_.coefficients(), localView_, localCoefficients_, Dune::PriorityTag<4>{});
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
98
99
100
101
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
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
    bound_ = true;
  }

  /// \brief Unbind the LocalView from the element
  void unbind()
  {
    localView_.unbind();
    bound_ = false;
  }

  /// \brief Evaluate LocalFunction at bound element in local coordinates
  Range operator()(Domain const& x) const
  {
    assert( bound_ );
    Range y(0);

    auto&& nodeToRangeEntry = globalFunction_.nodeToRangeEntry_;
    for_each_leaf_node(*subTree_, [&,this](auto const& node, auto const& tp)
    {
      auto localBasisCache = makeNodeCache(node);
      auto const& shapeFunctionValues = localBasisCache.evaluateFunction(localView_.element().type(), x);
      std::size_t size = node.finiteElement().size();

      // Get range entry associated to this node
      auto re = Dune::Functions::flatVectorView(nodeToRangeEntry(node, tp, y));

      for (std::size_t i = 0; i < size; ++i) {
        // Get coefficient associated to i-th shape function
        auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);

        // Get value of i-th shape function
        auto v = Dune::Functions::flatVectorView(shapeFunctionValues[i]);

        std::size_t dimC = c.size();
        std::size_t dimV = v.size();
        assert(dimC*dimV == std::size_t(re.size()));
        for(std::size_t j = 0; j < dimC; ++j) {
          auto&& c_j = c[j];
          for(std::size_t k = 0; k < dimV; ++k)
            re[j*dimV + k] += c_j*v[k];
        }
      }
    });

    return y;
  }

  /// \brief Create a LocalFunction representing the gradient. \relates GradientLocalFunction
  GradientLocalFunction makeDerivative(tag::gradient type) const
  {
    return GradientLocalFunction{globalFunction_, type};
  }

  DivergenceLocalFunction makeDerivative(tag::divergence type) const
  {
    return DivergenceLocalFunction{globalFunction_, type};
  }

  PartialLocalFunction makeDerivative(tag::partial type) const
  {
    return PartialLocalFunction{globalFunction_, type};
  }

  /// \brief The \ref polynomialDegree() of the LocalFunctions
  int order() const
  {
    assert( bound_ );
    return polynomialDegree(*subTree_);
  }

  /// \brief Return the bound element
  Element const& localContext() const
  {
    assert( bound_ );
    return localView_.element();
  }

private:
  DiscreteFunction globalFunction_;
  LocalView localView_;
  SubTree const* subTree_;

Praetorius, Simon's avatar
Praetorius, Simon committed
152
  std::vector<ValueType> localCoefficients_;
153
154
155
156
  bool bound_ = false;
};


Praetorius, Simon's avatar
Praetorius, Simon committed
157
template <class Coeff, class GB, class TP>
158
  template <class Type>
Praetorius, Simon's avatar
Praetorius, Simon committed
159
class DiscreteFunction<Coeff const,GB,TP>::DerivativeLocalFunctionBase
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
{
  using R = typename DiscreteFunction::Range;
  using D = typename DiscreteFunction::Domain;
  using RawSignature = typename Dune::Functions::SignatureTraits<R(D)>::RawSignature;
  using Traits = DerivativeTraits<RawSignature, Type>;

public:
  using Domain = typename EntitySet::LocalCoordinate;
  using Range = typename Traits::Range;

  enum { hasDerivative = false };

private:
  using LocalView = typename GlobalBasis::LocalView;
  using Element = typename EntitySet::Element;
  using Geometry = typename Element::Geometry;

public:
  /// Constructor. Stores a copy of the DiscreteFunction.
  DerivativeLocalFunctionBase(DiscreteFunction const& globalFunction, Type const& type)
    : globalFunction_(globalFunction)
    , type_(type)
Praetorius, Simon's avatar
Praetorius, Simon committed
182
    , localView_(globalFunction_.basis().localView())
183
    , subTree_(&Dune::TypeTree::child(localView_.tree(), globalFunction_.treePath()))
184
185
186
187
188
189
  {}

  /// Copy constructor.
  DerivativeLocalFunctionBase(DerivativeLocalFunctionBase const& other)
    : globalFunction_(other.globalFunction_)
    , type_(other.type_)
Praetorius, Simon's avatar
Praetorius, Simon committed
190
    , localView_(globalFunction_.basis().localView())
191
    , subTree_(&Dune::TypeTree::child(localView_.tree(), globalFunction_.treePath()))
192
193
194
195
196
197
198
  {}

  void bind(Element const& element)
  {
    localView_.bind(element);
    geometry_.emplace(element.geometry());

Praetorius, Simon's avatar
Praetorius, Simon committed
199
    Impl::gather(globalFunction_.coefficients(), localView_, localCoefficients_, Dune::PriorityTag<4>{});
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
    bound_ = true;
  }

  void unbind()
  {
    localView_.unbind();
    geometry_.reset();
    bound_ = false;
  }

  int order() const
  {
    assert( bound_ );
    return std::max(0, polynomialDegree(*subTree_)-1);
  }

  /// Return the bound element
  Element const& localContext() const
  {
    assert( bound_ );
    return localView_.element();
  }

  Geometry const& geometry() const
  {
    assert( bound_ );
    return geometry_.value();
  }

  LocalView const& localView() const
  {
    return localView_;
  }

protected:
  DiscreteFunction globalFunction_;
  Type type_;
  LocalView localView_;
  SubTree const* subTree_;
  Dune::Std::optional<Geometry> geometry_;
Praetorius, Simon's avatar
Praetorius, Simon committed
240
  std::vector<ValueType> localCoefficients_;
241
242
243
244
  bool bound_ = false;
};


Praetorius, Simon's avatar
Praetorius, Simon committed
245
246
template <class Coeff, class GB, class TP>
class DiscreteFunction<Coeff const,GB,TP>::GradientLocalFunction
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
    : public DerivativeLocalFunctionBase<tag::gradient>
{
  using Super = DerivativeLocalFunctionBase<tag::gradient>;
public:
  using Range = typename Super::Range;
  using Domain = typename Super::Domain;

  // adopt constructor from base class
  using DerivativeLocalFunctionBase<tag::gradient>::DerivativeLocalFunctionBase;

  /// Evaluate Gradient at bound element in local coordinates
  Range operator()(Domain const& x) const
  {
    assert( this->bound_ );
    Range dy(0);

    auto&& nodeToRangeEntry = this->globalFunction_.nodeToRangeEntry_;
    for_each_leaf_node(*this->subTree_, [&](auto const& node, auto const& tp)
    {
      auto localBasis = makeLocalToGlobalBasisAdapter(node, this->geometry());
      auto const& gradients = localBasis.gradientsAt(x);

      // Get range entry associated to this node
      auto re = Dune::Functions::flatVectorView(nodeToRangeEntry(node, tp, dy));

      for (std::size_t i = 0; i < localBasis.size(); ++i) {
        // Get coefficient associated to i-th shape function
        auto c = Dune::Functions::flatVectorView(this->localCoefficients_[node.localIndex(i)]);

        // Get value of i-th transformed reference gradient
        auto grad = Dune::Functions::flatVectorView(gradients[i]);

        std::size_t dimC = c.size();
        std::size_t dimV = grad.size();
        assert(dimC*dimV == std::size_t(re.size()));
        for(std::size_t j = 0; j < dimC; ++j) {
          auto&& c_j = c[j];
          for(std::size_t k = 0; k < dimV; ++k)
            re[j*dimV + k] += c_j*grad[k];
        }
      }
    });

    return dy;
  }

  using Super::localCoefficients_;
};


Praetorius, Simon's avatar
Praetorius, Simon committed
297
298
template <class Coeff, class GB, class TP>
class DiscreteFunction<Coeff const,GB,TP>::DivergenceLocalFunction
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
    : public DerivativeLocalFunctionBase<tag::divergence>
{
  using Super = DerivativeLocalFunctionBase<tag::divergence>;

public:
  using Range = typename Super::Range;
  using Domain = typename Super::Domain;

  // adopt constructor from base class
  using DerivativeLocalFunctionBase<tag::divergence>::DerivativeLocalFunctionBase;

  /// Evaluate divergence at bound element in local coordinates
  Range operator()(Domain const& x) const
  {
    return evaluate_node(x, bool_t<SubTree::isPower && SubTree::degree() == GridView::dimensionworld>{});
  }

private:
  Range evaluate_node(Domain const& x, std::false_type) const
  {
    error_exit("Cannot calculate divergence(node) where node is not a power node.");
    return Range(0);
  }

  Range evaluate_node(Domain const& x, std::true_type) const
  {
    static_assert(static_size_v<Range> == 1, "Implemented for scalar coefficients only.");

    assert( this->bound_ );
    Range dy(0);

    auto&& node = *this->subTree_;

    auto localBasis = makeLocalToGlobalBasisAdapter(node.child(0), this->geometry());
    auto const& gradients = localBasis.gradientsAt(x);

    auto re = Dune::Functions::flatVectorView(dy);
    assert(int(re.size()) == 1);
    for (std::size_t i = 0; i < localBasis.size(); ++i) {
      auto grad = Dune::Functions::flatVectorView(gradients[i]);

      assert(int(grad.size()) == GridView::dimensionworld);
      for (std::size_t j = 0; j < GridView::dimensionworld; ++j)
        re[0] += localCoefficients_[node.child(j).localIndex(i)] * grad[j];
    }

    return dy;
  }

  using Super::localCoefficients_;
};


Praetorius, Simon's avatar
Praetorius, Simon committed
352
353
template <class Coeff, class GB, class TP>
class DiscreteFunction<Coeff const,GB,TP>::PartialLocalFunction
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
    : public DerivativeLocalFunctionBase<tag::partial>
{
  using Super = DerivativeLocalFunctionBase<tag::partial>;

public:
  using Range = typename Super::Range;
  using Domain = typename Super::Domain;

  using DerivativeLocalFunctionBase<tag::partial>::DerivativeLocalFunctionBase;

  /// Evaluate partial derivative at bound element in local coordinates
  Range operator()(Domain const& x) const
  {
    assert( this->bound_ );
    Range dy(0);

    std::size_t comp = this->type_.comp;

    auto&& nodeToRangeEntry = this->globalFunction_.nodeToRangeEntry_;
    for_each_leaf_node(*this->subTree_, [&](auto const& node, auto const& tp)
    {
      auto localBasis = makeLocalToGlobalBasisAdapter(node, this->geometry());
      auto const& partial = localBasis.partialsAt(x, comp);

      // Get range entry associated to this node
      auto re = Dune::Functions::flatVectorView(nodeToRangeEntry(node, tp, dy));

      for (std::size_t i = 0; i < localBasis.size(); ++i) {
        // Get coefficient associated to i-th shape function
        auto c = Dune::Functions::flatVectorView(this->localCoefficients_[node.localIndex(i)]);

        // Get value of i-th transformed reference partial_derivative
        auto d_comp = Dune::Functions::flatVectorView(partial[i]);

        std::size_t dimC = c.size();
        std::size_t dimV = d_comp.size();
        assert(dimC*dimV == std::size_t(re.size()));
        for(std::size_t j = 0; j < dimC; ++j) {
          auto&& c_j = c[j];
          for(std::size_t k = 0; k < dimV; ++k)
            re[j*dimV + k] += c_j*d_comp[k];
        }
      }
    });

    return dy;
  }

  using Super::localCoefficients_;
};


} // end namespace AMDiS