DiscreteFunction.inc.hpp 12.7 KB
Newer Older
1
2
#pragma once

3
#include <amdis/common/DerivativeTraits.hpp>
4
#include <amdis/common/FieldMatVec.hpp>
5
#include <amdis/utility/LocalBasisCache.hpp>
6
#include <amdis/utility/LocalToGlobalAdapter.hpp>
7

8
#include <dune/common/ftraits.hh>
9
10
#include <dune/grid/utility/hierarchicsearch.hh>

11
12
namespace AMDiS {

13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
template <class GB, class VT, class TP>
class DiscreteFunction<GB,VT,TP>::LocalFunction
{
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:
28
  /// Constructor. Stores a copy of the DiscreteFunction.
29
  LocalFunction(DiscreteFunction const& globalFunction)
30
    : globalFunction_(globalFunction)
31
    , localView_(globalFunction_.basis()->localView())
32
    , subTree_(&child(localView_.tree(), globalFunction_.treePath()))
33
34
  {}

35
  /// Copy constructor.
36
37
  LocalFunction(LocalFunction const& other)
    : globalFunction_(other.globalFunction_)
38
    , localView_(globalFunction_.basis()->localView())
39
40
41
    , subTree_(&child(localView_.tree(), globalFunction_.treePath()))
  {}

42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
  /// \brief Bind the LocalView to the element
  void bind(Element const& element)
  {
    localView_.bind(element);
    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;

  /// \brief Create a LocalFunction representing the gradient. \relates GradientLocalFunction
60
  GradientLocalFunction makeDerivative(tag::gradient type) const
61
  {
62
63
64
65
66
67
68
69
70
71
72
    return GradientLocalFunction{globalFunction_, type};
  }

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

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

  /// \brief The \ref polynomialDegree() of the LocalFunctions
76
  int order() const
77
  {
78
79
    assert( bound_ );
    return polynomialDegree(*subTree_);
80
81
82
83
84
85
86
87
88
89
  }

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

private:
90
  DiscreteFunction globalFunction_;
91
92
93
94
95
96
97
  LocalView localView_;
  SubTree const* subTree_;
  bool bound_ = false;
};


template <class GB, class VT, class TP>
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
typename DiscreteFunction<GB,VT,TP>::Range DiscreteFunction<GB,VT,TP>::
LocalFunction::operator()(Domain const& x) const
{
  assert( bound_ );
  Range y(0);

  auto&& coefficients = *globalFunction_.dofVector_;
  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) {
      auto&& multiIndex = localView_.index(node.localIndex(i));

      // Get coefficient associated to i-th shape function
      auto c = Dune::Functions::flatVectorView(coefficients[multiIndex]);

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


template <class GB, class VT, class TP>
  template <class Type>
class DiscreteFunction<GB,VT,TP>::DerivativeLocalFunctionBase
142
143
144
145
{
  using R = typename DiscreteFunction::Range;
  using D = typename DiscreteFunction::Domain;
  using RawSignature = typename Dune::Functions::SignatureTraits<R(D)>::RawSignature;
146
  using Traits = DerivativeTraits<RawSignature, Type>;
147
148
149

public:
  using Domain = typename EntitySet::LocalCoordinate;
150
  using Range = typename Traits::Range;
151
152
153
154
155
156
157
158
159

  enum { hasDerivative = false };

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

public:
160
  /// Constructor. Stores a copy of the DiscreteFunction.
161
  DerivativeLocalFunctionBase(DiscreteFunction const& globalFunction, Type const& type)
162
    : globalFunction_(globalFunction)
163
    , type_(type)
164
    , localView_(globalFunction_.basis()->localView())
165
    , subTree_(&child(localView_.tree(), globalFunction_.treePath()))
166
167
  {}

168
  /// Copy constructor.
169
  DerivativeLocalFunctionBase(DerivativeLocalFunctionBase const& other)
170
    : globalFunction_(other.globalFunction_)
171
    , type_(other.type_)
172
    , localView_(globalFunction_.basis()->localView())
173
174
175
    , subTree_(&child(localView_.tree(), globalFunction_.treePath()))
  {}

176
177
178
179
180
181
182
183
184
185
186
187
188
189
  void bind(Element const& element)
  {
    localView_.bind(element);
    geometry_.emplace(element.geometry());
    bound_ = true;
  }

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

190
  int order() const
191
  {
192
193
    assert( bound_ );
    return std::max(0, polynomialDegree(*subTree_)-1);
194
195
196
197
198
199
200
201
202
  }

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

203
protected:
204
  DiscreteFunction globalFunction_;
205
  Type type_;
206
207
208
209
210
211
212
213
  LocalView localView_;
  SubTree const* subTree_;
  Dune::Std::optional<Geometry> geometry_;
  bool bound_ = false;
};


template <class GB, class VT, class TP>
214
215
class DiscreteFunction<GB,VT,TP>::GradientLocalFunction
    : public DerivativeLocalFunctionBase<tag::gradient>
216
{
217
218
219
220
  using Super = DerivativeLocalFunctionBase<tag::gradient>;
public:
  using Range = typename Super::Range;
  using Domain = typename Super::Domain;
221

222
223
224
225
  using Super::DerivativeLocalFunctionBase;

  /// Evaluate Gradient at bound element in local coordinates
  Range operator()(Domain const& x) const
226
  {
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
    assert( Super::bound_ );
    Range dy(0);

    auto&& coefficients = *Super::globalFunction_.dofVector_;
    auto&& nodeToRangeEntry = Super::globalFunction_.nodeToRangeEntry_;
    for_each_leaf_node(*Super::subTree_, [&,this](auto const& node, auto const& tp)
    {
      auto localBasis = makeLocalToGlobalBasisAdapter(node, Super::geometry_.value());
      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) {
        auto&& multiIndex = Super::localView_.index(node.localIndex(i));

        // Get coefficient associated to i-th shape function
        auto c = Dune::Functions::flatVectorView(coefficients[multiIndex]);

        // 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];
        }
      }
    });
259

260
261
262
    return dy;
  }
};
263
264


265
266
267
268
269
template <class GB, class VT, class TP>
class DiscreteFunction<GB,VT,TP>::DivergenceLocalFunction
    : public DerivativeLocalFunctionBase<tag::divergence>
{
  using Super = DerivativeLocalFunctionBase<tag::divergence>;
270

271
272
273
public:
  using Range = typename Super::Range;
  using Domain = typename Super::Domain;
274

275
  using Super::DerivativeLocalFunctionBase;
276

277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
  /// 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(Size_v<Range> == 1, "Implemented for scalar coefficients only.");

    assert( Super::bound_ );
    Range dy(0);

    auto&& coefficients = *Super::globalFunction_.dofVector_;
    auto&& node = *Super::subTree_;

    auto localBasis = makeLocalToGlobalBasisAdapter(node.child(0), Super::geometry_.value());
    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) {
        auto&& multiIndex = Super::localView_.index(node.child(j).localIndex(i));
        re[0] += coefficients[multiIndex] * grad[j];
312
313
314
      }
    }

315
316
317
    return dy;
  }
};
318

Praetorius, Simon's avatar
Praetorius, Simon committed
319

320
template <class GB, class VT, class TP>
321
322
class DiscreteFunction<GB,VT,TP>::PartialLocalFunction
    : public DerivativeLocalFunctionBase<tag::partial>
323
{
324
  using Super = DerivativeLocalFunctionBase<tag::partial>;
325

326
327
328
public:
  using Range = typename Super::Range;
  using Domain = typename Super::Domain;
329

330
  using Super::DerivativeLocalFunctionBase;
331

332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
  /// Evaluate partial derivative at bound element in local coordinates
  Range operator()(Domain const& x) const
  {
    assert( Super::bound_ );
    Range dy(0);

    auto&& coefficients = *Super::globalFunction_.dofVector_;
    auto&& nodeToRangeEntry = Super::globalFunction_.nodeToRangeEntry_;
    for_each_leaf_node(*Super::subTree_, [&,this](auto const& node, auto const& tp)
    {
      auto localBasis = makeLocalToGlobalBasisAdapter(node, Super::geometry_.value());
      auto const& partial = localBasis.partialsAt(x, Super::type_.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) {
        auto&& multiIndex = Super::localView_.index(node.localIndex(i));

        // Get coefficient associated to i-th shape function
        auto c = Dune::Functions::flatVectorView(coefficients[multiIndex]);

        // 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];
        }
      }
    });
367

368
369
    return dy;
  }
370

371
372
373
374
375
376
377
private:
  template <class T, class J, class RG>
  void test_partial(std::vector<T> const& partial, std::size_t comp, J const& jacobian, RG const& referenceGradients) const
  {
    for (std::size_t i = 0; i < referenceGradients.size(); ++i) { // J^(-T) * D[phi] -> grad^T
      auto grad = jacobian * Dune::MatVec::as_vector(referenceGradients[i]);
      auto grad_ = Dune::Functions::flatVectorView(grad);
378

379
380
381
382
383
384
      if (std::abs(grad_[comp] - partial[i]) > 1.e-10) {
        print_error(partial[i], comp, grad, jacobian, referenceGradients[i][0]);
        error_exit("wrong partial derivatives");
      }
    }
  }
385

386
387
388
389
390
391
392
393
394
395
  template <class T, class G, class J, class RG>
  void print_error(T const& partial_i, std::size_t comp, G const& grad_i, J const& jacobian, RG const& referenceGradient) const
  {
    msg(__PRETTY_FUNCTION__);
    msg("jacobian = {}", jacobian);
    auto&& jacobian_mat = Dune::MatVec::as_matrix(jacobian);
    msg("jacobian_ = [");
    for (std::size_t i = 0; i < jacobian_mat.N(); ++i) {
      for (std::size_t j = 0; j < jacobian_mat.M(); ++j) {
        msg_("{} ",jacobian_mat[i][j]);
396
      }
397
      msg(";");
398
    }
399
    msg("]");
400

401
402
403
404
405
    msg("ref-grad = {}", referenceGradient);
    auto grad_ = Dune::Functions::flatVectorView(grad_i);
    msg("comp = {}, grad_ = ({}, {}), partial = {}", comp, grad_[0], grad_[1], partial_i);
  }
};
406

407

408
409
template <class GB, class VT, class TP>
typename DiscreteFunction<GB,VT,TP>::Range DiscreteFunction<GB,VT,TP>::
410
411
operator()(Domain const& x) const
{
412
413
  using Grid = typename GlobalBasis::GridView::Grid;
  using IS = typename GlobalBasis::GridView::IndexSet;
414

415
  auto const& gv = this->basis()->gridView();
416
417
418
419
420
421
422
423
424
  Dune::HierarchicSearch<Grid,IS> hsearch{gv.grid(), gv.indexSet()};

  auto element = hsearch.findEntity(x);
  auto geometry = element.geometry();
  auto localFct = localFunction(*this);
  localFct.bind(element);
  return localFct(geometry.local(x));
}

425
} // end namespace AMDiS