lagrangegridcreator.hh 14.2 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
3
4
5
#pragma once

#include <cassert>
#include <cstdint>
#include <limits>
6
#include <optional>
Praetorius, Simon's avatar
Praetorius, Simon committed
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
#include <vector>

#include <dune/common/exceptions.hh>
#include <dune/common/hybridutilities.hh>
#include <dune/geometry/utility/typefromvertexcount.hh>
#include <dune/localfunctions/lagrange.hh>
#include <dune/grid/common/gridfactory.hh>

#include <dune/vtk/forward.hh>
#include <dune/vtk/vtktypes.hh>
#include <dune/vtk/gridcreatorinterface.hh>
#include <dune/vtk/utility/lagrangepoints.hh>

namespace Dune
{
  // \brief Create a grid from data that represents higher (lagrange) cells.
  /**
24
   * The grid is created from the first nodes of a cell parametrization, representing
Praetorius, Simon's avatar
Praetorius, Simon committed
25
26
27
28
29
30
   * the  corner vertices. Thus a piecewise "flat" grid is constructed. The
   * parametrization is 1. passed as a local element parametrization to the
   * `insertElement()` function of a gridFactory to allow the grid itself to handle the
   * parametrization and 2. is stored internally that can be accessed by using this
   * GridCreator object as a grid function, or by extracting locally the parametrization
   * on each existing grid element after creation of the grid.
31
   *
Praetorius, Simon's avatar
Praetorius, Simon committed
32
33
   * So, the LagrangeGridCreator models both, a `GridCreator` and a `GridFunction`.
   **/
Praetorius, Simon's avatar
Praetorius, Simon committed
34
  template <class GridType>
Praetorius, Simon's avatar
Praetorius, Simon committed
35
  struct LagrangeGridCreator
Praetorius, Simon's avatar
Praetorius, Simon committed
36
      : public GridCreatorInterface<GridType, LagrangeGridCreator<GridType>>
Praetorius, Simon's avatar
Praetorius, Simon committed
37
38
  {
    using Self = LagrangeGridCreator;
Praetorius, Simon's avatar
Praetorius, Simon committed
39
    using Super = GridCreatorInterface<GridType, Self>;
Praetorius, Simon's avatar
Praetorius, Simon committed
40
41
42
43
44
45
    using GlobalCoordinate = typename Super::GlobalCoordinate;

    using Nodes = std::vector<GlobalCoordinate>;

    struct ElementParametrization
    {
46
      GeometryType type;                  //< Geometry type of the element
Praetorius, Simon's avatar
Praetorius, Simon committed
47
48
49
50
51
      std::vector<std::int64_t> nodes;    //< Indices of the w.r.t. `nodes_` vector
      std::vector<unsigned int> corners;  //< Insertion-indices of the element corner nodes
    };

    using Parametrization = std::vector<ElementParametrization>;
Praetorius, Simon's avatar
Praetorius, Simon committed
52
    using Element = typename GridType::template Codim<0>::Entity;
Praetorius, Simon's avatar
Praetorius, Simon committed
53
54
55
56
57
58
59
60
    using LocalCoordinate = typename Element::Geometry::LocalCoordinate;

    class LocalParametrization;
    class LocalFunction;

  public:
    using Super::factory;

Praetorius, Simon's avatar
Praetorius, Simon committed
61
    LagrangeGridCreator (GridFactory<GridType>& factory)
Praetorius, Simon's avatar
Praetorius, Simon committed
62
63
64
65
66
67
68
69
70
71
72
73
      : Super(factory)
    {}

    /// Implementation of the interface function `insertVertices()`
    void insertVerticesImpl (std::vector<GlobalCoordinate> const& points,
                             std::vector<std::uint64_t> const& /*point_ids*/)
    {
      // store point coordinates in member variable
      nodes_ = points;
    }

    template <class F>
74
    using HasParametrizedElements = decltype(std::declval<F>().insertElement(std::declval<GeometryType>(),
Praetorius, Simon's avatar
Praetorius, Simon committed
75
76
77
78
79
80
81
82
83
84
85
      std::declval<std::vector<unsigned int> const&>(), std::declval<std::function<GlobalCoordinate(LocalCoordinate)>>()));

    /// Implementation of the interface function `insertElements()`
    void insertElementsImpl (std::vector<std::uint8_t> const& types,
                             std::vector<std::int64_t> const& offsets,
                             std::vector<std::int64_t> const& connectivity)
    {
      assert(nodes_.size() > 0);

      // mapping of node index to element-vertex index
      std::vector<std::int64_t> elementVertices(nodes_.size(), -1);
86
      parametrization_.reserve(types.size());
Praetorius, Simon's avatar
Praetorius, Simon committed
87
88
89
90

      std::int64_t vertexIndex = 0;
      for (std::size_t i = 0; i < types.size(); ++i) {
        auto type = Vtk::to_geometry(types[i]);
Praetorius, Simon's avatar
Praetorius, Simon committed
91
        if (type.dim() != GridType::dimension)
92
93
          continue;

Praetorius, Simon's avatar
Praetorius, Simon committed
94
        Vtk::CellType cellType{type};
Praetorius, Simon's avatar
Praetorius, Simon committed
95
        auto refElem = referenceElement<double,GridType::dimension>(type);
Praetorius, Simon's avatar
Praetorius, Simon committed
96
97
98

        std::int64_t shift = (i == 0 ? 0 : offsets[i-1]);
        int nNodes = offsets[i] - shift;
Praetorius, Simon's avatar
Praetorius, Simon committed
99
        int nVertices = refElem.size(GridType::dimension);
Praetorius, Simon's avatar
Praetorius, Simon committed
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122

        // insert vertices into grid and construct element vertices
        std::vector<unsigned int> element(nVertices);
        for (int j = 0; j < nVertices; ++j) {
          auto index = connectivity.at(shift + j);
          auto& vertex = elementVertices.at(index);
          if (vertex < 0) {
            factory().insertVertex(nodes_.at(index));
            vertex = vertexIndex++;
          }
          element[j] = vertex;
        }

        // permute element indices
        if (!cellType.noPermutation()) {
          // apply index permutation
          std::vector<unsigned int> cell(element.size());
          for (int j = 0; j < element.size(); ++j)
            cell[j] = element[cellType.permutation(j)];
          std::swap(element, cell);
        }

        // fill vector of element parametrizations
123
124
125
126
        parametrization_.push_back(ElementParametrization{type});
        auto& param = parametrization_.back();

        param.nodes.resize(nNodes);
Praetorius, Simon's avatar
Praetorius, Simon committed
127
        for (int j = 0; j < nNodes; ++j)
128
129
          param.nodes[j] = connectivity.at(shift + j);
        param.corners = element;
Praetorius, Simon's avatar
Praetorius, Simon committed
130
131

        // try to create element with parametrization
Praetorius, Simon's avatar
Praetorius, Simon committed
132
        if constexpr (Std::is_detected_v<HasParametrizedElements, GridFactory<GridType>>) {
Praetorius, Simon's avatar
Praetorius, Simon committed
133
          try {
134
135
            factory().insertElement(type, element,
              localParametrization(parametrization_.size()-1));
Praetorius, Simon's avatar
Praetorius, Simon committed
136
137
138
139
140
141
142
143
144
145
146
147
          } catch (Dune::GridError const& /* notImplemented */) {
            factory().insertElement(type, element);
          }
        } else {
          factory().insertElement(type, element);
        }
      }
    }

    /// \brief Construct an element parametrization
    /**
     * The returned LocalParametrization is a mapping `GlobalCoordinate(LocalCoordinate)`
148
149
     * where `LocalCoordinate is w.r.t. the local coordinate system in an element with
     * given `insertionIndex` (defined by the inserted corner vertices) and
Praetorius, Simon's avatar
Praetorius, Simon committed
150
151
     * `GlobalCoordinate` a world coordinate in the parametrized grid.
     **/
152
    LocalParametrization localParametrization (unsigned int insertionIndex) const
Praetorius, Simon's avatar
Praetorius, Simon committed
153
154
155
156
157
158
159
160
161
    {
      assert(!nodes_.empty() && !parametrization_.empty());
      auto const& localParam = parametrization_.at(insertionIndex);
      return LocalParametrization{nodes_, localParam, order(localParam)};
    }

    /// \brief Construct an element parametrization
    /**
     * The returned LocalParametrization is a mapping `GlobalCoordinate(LocalCoordinate)`
162
     * where `LocalCoordinate is w.r.t. the local coordinate system in the passed element
Praetorius, Simon's avatar
Praetorius, Simon committed
163
     * and `GlobalCoordinate` a world coordinate in the parametrized grid.
164
     *
Praetorius, Simon's avatar
Praetorius, Simon committed
165
     * Note, when an element is passed, it might have a different local coordinate system
Praetorius, Simon's avatar
Praetorius, Simon committed
166
167
168
169
170
     * than the coordinate system used to defined the element parametrization. Thus
     * coordinate transform is internally chained to the evaluation of the local
     * parametrization. This local geometry transform is obtained by figuring out the
     * permutation of corners in the element corresponding to the inserted corner
     * vertices.
Praetorius, Simon's avatar
Praetorius, Simon committed
171
172
173
174
175
176
177
178
179
180
     **/
    LocalParametrization localParametrization (Element const& element) const
    {
      assert(!nodes_.empty() && !parametrization_.empty());

      unsigned int insertionIndex = factory().insertionIndex(element);
      auto const& localParam = parametrization_.at(insertionIndex);
      assert(element.type() == localParam.type);

      // collect indices of vertices
Praetorius, Simon's avatar
Praetorius, Simon committed
181
182
183
      std::vector<unsigned int> indices(element.subEntities(GridType::dimension));
      for (unsigned int i = 0; i < element.subEntities(GridType::dimension); ++i)
        indices[i] = factory().insertionIndex(element.template subEntity<GridType::dimension>(i));
Praetorius, Simon's avatar
Praetorius, Simon committed
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209

      // calculate permutation vector
      std::vector<unsigned int> permutation(indices.size());
      for (std::size_t i = 0; i < indices.size(); ++i) {
        auto it = std::find(localParam.corners.begin(), localParam.corners.end(), indices[i]);
        assert(it != localParam.corners.end());
        permutation[i] = std::distance(localParam.corners.begin(), it);
      }

      return LocalParametrization{nodes_, localParam, order(localParam), permutation};
    }

    /// Determine lagrange order from number of points
    template <class LocalParam>
    unsigned int order (LocalParam const& localParam) const
    {
      GeometryType type = localParam.type;
      std::size_t nNodes = localParam.nodes.size();
      for (unsigned int o = 1; o <= nNodes; ++o)
        if (numLagrangePoints(type.id(), type.dim(), o) == nNodes)
          return o;

      return 1;
    }

    /// Determine lagrange order from number of points from the first element parametrization
210
    unsigned int order () const
Praetorius, Simon's avatar
Praetorius, Simon committed
211
212
213
214
215
    {
      assert(!parametrization_.empty());
      return order(parametrization_.front());
    }

216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
  public:
    /// \brief Local function representing the parametrization of the grid.
    /**
     * The returned object models Functions::Concept::LocalFunction
     * and can thus be bound to an element of the created grid and evaluated in
     * the local coordinates of the bound element.
     *
     * It is implemented in terms of the \ref LocalParametrization function
     * returned by the method \ref localParametrization(element). See comments
     * there for further details.
     *
     * Note, this methods requires the GridCreator to be based by
     * lvalue-reference. This is necessary, since we want to guarantee that all
     * internal storage is preserved while evaluating the local function.
     **/
    friend LocalFunction localFunction (LagrangeGridCreator& gridCreator)
    {
      return LocalFunction{gridCreator};
    }

236
237
238
239
240
241
242
243
244
245
246
    friend LocalFunction localFunction (LagrangeGridCreator const& gridCreator)
    {
      return LocalFunction{gridCreator};
    }

    friend LocalFunction localFunction (LagrangeGridCreator&& gridCreator)
    {
      DUNE_THROW(Dune::Exception, "Cannot pass temporary LagrangeGridCreator to localFunction(). Pass an lvalue-reference instead.");
      return LocalFunction{gridCreator};
    }

247
248
    struct EntitySet
    {
Praetorius, Simon's avatar
Praetorius, Simon committed
249
      using Grid = GridType;
250
      using GlobalCoordinate = typename Self::GlobalCoordinate;
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
    };

    /// Dummy function returning a placeholder entityset
    EntitySet entitySet () const
    {
      assert(false && "Should not be used!");
      return EntitySet{};
    }

    /// Dummy function returning a placeholder entityset
    GlobalCoordinate operator() (GlobalCoordinate const&) const
    {
      assert(false && "Should not be used!");
      return GlobalCoordinate{};
    }

Praetorius, Simon's avatar
Praetorius, Simon committed
267
268
269
270
271
272
273
274
275
276
277
278
279
  private:
    /// All point coordinates inclusing the higher-order lagrange points
    Nodes nodes_;

    /// Parametrization for all elements
    Parametrization parametrization_;
  };


  template <class Grid>
  class LagrangeGridCreator<Grid>::LocalParametrization
  {
    using ctype = typename Grid::ctype;
280

Praetorius, Simon's avatar
Praetorius, Simon committed
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
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
    using GlobalCoordinate = typename Grid::template Codim<0>::Entity::Geometry::GlobalCoordinate;
    using LocalCoordinate = typename Grid::template Codim<0>::Entity::Geometry::LocalCoordinate;
    using LocalGeometry = MultiLinearGeometry<ctype,Grid::dimension,Grid::dimension>;

    using LocalFE = LagrangeLocalFiniteElement<VtkLagrangePointSet, Grid::dimension, ctype, ctype>;
    using LocalBasis = typename LocalFE::Traits::LocalBasisType;
    using LocalBasisTraits = typename LocalBasis::Traits;

  public:
    /// Construct a local element parametrization
    template <class Nodes, class LocalParam>
    LocalParametrization (Nodes const& nodes, LocalParam const& param, unsigned int order)
      : localFE_(param.type, order)
      , localNodes_(param.nodes.size())
    {
      for (std::size_t i = 0; i < localNodes_.size(); ++i)
        localNodes_[i] = nodes[param.nodes[i]];
    }

    /// Construct a local element parametrization for elements with permuted corners
    template <class Nodes, class LocalParam, class Permutation>
    LocalParametrization (Nodes const& nodes, LocalParam const& param, unsigned int order, Permutation const& permutation)
      : LocalParametrization(nodes, param, order)
    {
      auto refElem = referenceElement<ctype,Grid::dimension>(param.type);
      std::vector<LocalCoordinate> corners(permutation.size());
      for (std::size_t i = 0; i < permutation.size(); ++i)
        corners[i] = refElem.position(permutation[i], Grid::dimension);

      localGeometry_.emplace(param.type, corners);
    }

    /// Evaluate the local parametrization in local coordinates
    template <class LocalCoordinate>
    GlobalCoordinate operator() (LocalCoordinate const& local) const
    {
      // map coordinates if element corners are permuted
      LocalCoordinate x = localGeometry_ ? localGeometry_->global(local) : local;

      LocalBasis const& localBasis = localFE_.localBasis();
      localBasis.evaluateFunction(x, shapeValues_);
      assert(shapeValues_.size() == localNodes_.size());

      GlobalCoordinate out(0);
      for (std::size_t i = 0; i < shapeValues_.size(); ++i)
        out.axpy(shapeValues_[i], localNodes_[i]);

      return out;
    }

  private:
    LocalFE localFE_;
    std::vector<GlobalCoordinate> localNodes_;
334
    std::optional<LocalGeometry> localGeometry_;
Praetorius, Simon's avatar
Praetorius, Simon committed
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349

    mutable std::vector<typename LocalBasisTraits::RangeType> shapeValues_;
  };


  template <class Grid>
  class LagrangeGridCreator<Grid>::LocalFunction
  {
    using ctype = typename Grid::ctype;
    using LocalContext = typename Grid::template Codim<0>::Entity;
    using GlobalCoordinate = typename LocalContext::Geometry::GlobalCoordinate;
    using LocalCoordinate = typename LocalContext::Geometry::LocalCoordinate;
    using LocalParametrization = typename LagrangeGridCreator::LocalParametrization;

  public:
350
    explicit LocalFunction (LagrangeGridCreator const& gridCreator)
Praetorius, Simon's avatar
Praetorius, Simon committed
351
352
353
      : gridCreator_(&gridCreator)
    {}

354
355
    explicit LocalFunction (LagrangeGridCreator&& gridCreator) = delete;

Praetorius, Simon's avatar
Praetorius, Simon committed
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
    /// Collect a local parametrization on the element
    void bind (LocalContext const& element)
    {
      localContext_ = element;
      localParametrization_.emplace(gridCreator_->localParametrization(element));
    }

    void unbind () { /* do nothing */ }

    /// Evaluate the local parametrization in local coordinates
    GlobalCoordinate operator() (LocalCoordinate const& local) const
    {
      assert(!!localParametrization_);
      return (*localParametrization_)(local);
    }

    /// Return the bound element
    LocalContext const& localContext () const
    {
      return localContext_;
    }

  private:
    LagrangeGridCreator const* gridCreator_;

    LocalContext localContext_;
382
    std::optional<LocalParametrization> localParametrization_;
Praetorius, Simon's avatar
Praetorius, Simon committed
383
384
385
  };

} // end namespace Dune