Dune Core Modules (2.9.0)

lagrangebasis.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_LAGRANGEBASIS_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_LAGRANGEBASIS_HH
5
6#include <type_traits>
8
10#include <dune/localfunctions/lagrange/equidistantpoints.hh>
11#include <dune/localfunctions/lagrange/pqkfactory.hh>
12
13#include <dune/functions/functionspacebases/nodes.hh>
14#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
15
16
17namespace Dune {
18namespace Functions {
19
20// *****************************************************************************
21// This is the reusable part of the LagrangeBasis. It contains
22//
23// LagrangePreBasis
24// LagrangeNode
25//
26// The pre-basis allows to create the others and is the owner of possible shared
27// state. These components do _not_ depend on the global basis and local view
28// and can be used without a global basis.
29// *****************************************************************************
30
31template<typename GV, int k, typename R=double>
32class LagrangeNode;
33
34template<typename GV, int k, typename R=double>
35class LagrangePreBasis;
36
37
38
53template<typename GV, int k, typename R>
55{
56 static const int dim = GV::dimension;
57 static const bool useDynamicOrder = (k<0);
58
59public:
60
62 using GridView = GV;
63
65 using size_type = std::size_t;
66
68 using Node = LagrangeNode<GV, k, R>;
69
70 static constexpr size_type maxMultiIndexSize = 1;
71 static constexpr size_type minMultiIndexSize = 1;
72 static constexpr size_type multiIndexBufferSize = 1;
73
76 : LagrangePreBasis(gv, std::numeric_limits<unsigned int>::max())
77 {}
78
80 LagrangePreBasis(const GridView& gv, unsigned int order) :
81 gridView_(gv), order_(order)
82 {
83 if (!useDynamicOrder && order!=std::numeric_limits<unsigned int>::max())
84 DUNE_THROW(RangeError, "Template argument k has to be -1 when supplying a run-time order!");
85
86 for (int i=0; i<=dim; i++)
87 {
88 dofsPerCube_[i] = computeDofsPerCube(i);
89 dofsPerSimplex_[i] = computeDofsPerSimplex(i);
90 }
91 dofsPerPrism_ = computeDofsPerPrism();
92 dofsPerPyramid_ = computeDofsPerPyramid();
93 }
94
97 {
98 vertexOffset_ = 0;
99 edgeOffset_ = vertexOffset_ + dofsPerCube(0) * ((size_type)gridView_.size(dim));
100
101 if (dim>=2)
102 {
103 triangleOffset_ = edgeOffset_ + dofsPerCube(1) * ((size_type) gridView_.size(dim-1));
104
105 quadrilateralOffset_ = triangleOffset_ + dofsPerSimplex(2) * ((size_type)gridView_.size(Dune::GeometryTypes::triangle));
106 }
107
108 if (dim==3) {
109 tetrahedronOffset_ = quadrilateralOffset_ + dofsPerCube(2) * ((size_type)gridView_.size(Dune::GeometryTypes::quadrilateral));
110
111 prismOffset_ = tetrahedronOffset_ + dofsPerSimplex(3) * ((size_type)gridView_.size(Dune::GeometryTypes::tetrahedron));
112
113 hexahedronOffset_ = prismOffset_ + dofsPerPrism() * ((size_type)gridView_.size(Dune::GeometryTypes::prism));
114
115 pyramidOffset_ = hexahedronOffset_ + dofsPerCube(3) * ((size_type)gridView_.size(Dune::GeometryTypes::hexahedron));
116 }
117 }
118
120 const GridView& gridView() const
121 {
122 return gridView_;
123 }
124
126 void update (const GridView& gv)
127 {
128 gridView_ = gv;
129 }
130
135 {
136 return Node{order_};
137 }
138
141 {
142 switch (dim)
143 {
144 case 1:
145 return dofsPerCube(0) * ((size_type)gridView_.size(dim))
146 + dofsPerCube(1) * ((size_type)gridView_.size(dim-1));
147 case 2:
148 {
149 return dofsPerCube(0) * ((size_type)gridView_.size(dim))
150 + dofsPerCube(1) * ((size_type)gridView_.size(dim-1))
153 }
154 case 3:
155 {
156 return dofsPerCube(0) * ((size_type)gridView_.size(dim))
157 + dofsPerCube(1) * ((size_type)gridView_.size(dim-1))
161 + dofsPerPyramid() * ((size_type)gridView_.size(Dune::GeometryTypes::pyramid))
162 + dofsPerPrism() * ((size_type)gridView_.size(Dune::GeometryTypes::prism))
163 + dofsPerCube(3) * ((size_type)gridView_.size(Dune::GeometryTypes::hexahedron));
164 }
165 }
166 DUNE_THROW(Dune::NotImplemented, "No size method for " << dim << "d grids available yet!");
167 }
168
170 template<class SizePrefix>
171 size_type size(const SizePrefix& prefix) const
172 {
173 assert(prefix.size() == 0 || prefix.size() == 1);
174 return (prefix.size() == 0) ? size() : 0;
175 }
176
179 {
180 return size();
181 }
182
185 {
186 // That cast to unsigned int is necessary because GV::dimension is an enum,
187 // which is not recognized by the power method as an integer type...
188 return power(order()+1, (unsigned int)GV::dimension);
189 }
190
191 template<typename It>
192 It indices(const Node& node, It it) const
193 {
194 for (size_type i = 0, end = node.finiteElement().size() ; i < end ; ++it, ++i)
195 {
196 Dune::LocalKey localKey = node.finiteElement().localCoefficients().localKey(i);
197 const auto& gridIndexSet = gridView().indexSet();
198 const auto& element = node.element();
199
200 // The dimension of the entity that the current dof is related to
201 auto dofDim = dim - localKey.codim();
202
203 // Test for a vertex dof
204 // The test for k==1 is redundant, but having it here allows the compiler to conclude
205 // at compile-time that the dofDim==0 case is the only one that will ever happen.
206 // This leads to measurable speed-up: see
207 // https://gitlab.dune-project.org/staging/dune-functions/issues/30
208 if (k==1 || dofDim==0) {
209 *it = {{ (size_type)(gridIndexSet.subIndex(element,localKey.subEntity(),dim)) }};
210 continue;
211 }
212
213 if (dofDim==1)
214 { // edge dof
215 if (dim==1) // element dof -- any local numbering is fine
216 {
217 *it = {{ edgeOffset_
218 + dofsPerCube(1) * ((size_type)gridIndexSet.subIndex(element,0,0))
219 + localKey.index() }};
220 continue;
221 }
222 else
223 {
224 const auto refElement
225 = Dune::referenceElement<double,dim>(element.type());
226
227 // We have to reverse the numbering if the local element edge is
228 // not aligned with the global edge.
229 auto v0 = (size_type)gridIndexSet.subIndex(element,refElement.subEntity(localKey.subEntity(),localKey.codim(),0,dim),dim);
230 auto v1 = (size_type)gridIndexSet.subIndex(element,refElement.subEntity(localKey.subEntity(),localKey.codim(),1,dim),dim);
231 bool flip = (v0 > v1);
232 *it = {{ (flip)
233 ? edgeOffset_
234 + dofsPerCube(1)*((size_type)gridIndexSet.subIndex(element,localKey.subEntity(),localKey.codim()))
235 + (dofsPerCube(1)-1)-localKey.index()
236 : edgeOffset_
237 + dofsPerCube(1)*((size_type)gridIndexSet.subIndex(element,localKey.subEntity(),localKey.codim()))
238 + localKey.index() }};
239 continue;
240 }
241 }
242
243 if (dofDim==2)
244 {
245 if (dim==2) // element dof -- any local numbering is fine
246 {
247 if (element.type().isTriangle())
248 {
249 *it = {{ triangleOffset_ + dofsPerSimplex(2)*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
250 continue;
251 }
252 else if (element.type().isQuadrilateral())
253 {
254 *it = {{ quadrilateralOffset_ + dofsPerCube(2)*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
255 continue;
256 }
257 else
258 DUNE_THROW(Dune::NotImplemented, "2d elements have to be triangles or quadrilaterals");
259 } else
260 {
261 const auto refElement
262 = Dune::referenceElement<double,dim>(element.type());
263
264 if (order()>3)
265 DUNE_THROW(Dune::NotImplemented, "LagrangeBasis for 3D grids is only implemented if k<=3");
266
267 if (order()==3 and !refElement.type(localKey.subEntity(), localKey.codim()).isTriangle())
268 DUNE_THROW(Dune::NotImplemented, "LagrangeBasis for 3D grids with k==3 is only implemented if the grid is a simplex grid");
269
270 *it = {{ triangleOffset_ + ((size_type)gridIndexSet.subIndex(element,localKey.subEntity(),localKey.codim())) }};
271 continue;
272 }
273 }
274
275 if (dofDim==3)
276 {
277 if (dim==3) // element dof -- any local numbering is fine
278 {
279 if (element.type().isTetrahedron())
280 {
281 *it = {{ tetrahedronOffset_ + dofsPerSimplex(3)*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
282 continue;
283 }
284 else if (element.type().isHexahedron())
285 {
286 *it = {{ hexahedronOffset_ + dofsPerCube(3)*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
287 continue;
288 }
289 else if (element.type().isPrism())
290 {
291 *it = {{ prismOffset_ + dofsPerPrism()*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
292 continue;
293 }
294 else if (element.type().isPyramid())
295 {
296 *it = {{ pyramidOffset_ + dofsPerPyramid()*((size_type)gridIndexSet.subIndex(element,0,0)) + localKey.index() }};
297 continue;
298 }
299 else
300 DUNE_THROW(Dune::NotImplemented, "3d elements have to be tetrahedra, hexahedra, prisms, or pyramids");
301 } else
302 DUNE_THROW(Dune::NotImplemented, "Grids of dimension larger than 3 are no supported");
303 }
304 DUNE_THROW(Dune::NotImplemented, "Grid contains elements not supported for the LagrangeBasis");
305 }
306 return it;
307 }
308
310 unsigned int order() const
311 {
312 return (useDynamicOrder) ? order_ : k;
313 }
314
315protected:
316 GridView gridView_;
317
318 // Run-time order, only valid if k<0
319 const unsigned int order_;
320
322 size_type dofsPerSimplex(std::size_t simplexDim) const
323 {
324 return useDynamicOrder ? dofsPerSimplex_[simplexDim] : computeDofsPerSimplex(simplexDim);
325 }
326
328 size_type dofsPerCube(std::size_t cubeDim) const
329 {
330 return useDynamicOrder ? dofsPerCube_[cubeDim] : computeDofsPerCube(cubeDim);
331 }
332
333 size_type dofsPerPrism() const
334 {
335 return useDynamicOrder ? dofsPerPrism_ : computeDofsPerPrism();
336 }
337
338 size_type dofsPerPyramid() const
339 {
340 return useDynamicOrder ? dofsPerPyramid_ : computeDofsPerPyramid();
341 }
342
344 size_type computeDofsPerSimplex(std::size_t simplexDim) const
345 {
346 return order() == 0 ? (dim == simplexDim ? 1 : 0) : Dune::binomial(std::size_t(order()-1),simplexDim);
347 }
348
350 size_type computeDofsPerCube(std::size_t cubeDim) const
351 {
352 return order() == 0 ? (dim == cubeDim ? 1 : 0) : Dune::power(order()-1, cubeDim);
353 }
354
355 size_type computeDofsPerPrism() const
356 {
357 return order() == 0 ? (dim == 3 ? 1 : 0) : (order()-1)*(order()-1)*(order()-2)/2;
358 }
359
360 size_type computeDofsPerPyramid() const
361 {
362 return order() == 0 ? (dim == 3 ? 1 : 0) : (order()-2)*(order()-1)*(2*order()-3)/6;
363 }
364
365 // When the order is given at run-time, the following numbers are pre-computed:
366 std::array<size_type,dim+1> dofsPerSimplex_;
367 std::array<size_type,dim+1> dofsPerCube_;
368 size_type dofsPerPrism_;
369 size_type dofsPerPyramid_;
370
371 size_type vertexOffset_;
372 size_type edgeOffset_;
373 size_type triangleOffset_;
374 size_type quadrilateralOffset_;
375 size_type tetrahedronOffset_;
376 size_type pyramidOffset_;
377 size_type prismOffset_;
378 size_type hexahedronOffset_;
379
380};
381
382
383
384template<typename GV, int k, typename R>
385class LagrangeNode :
386 public LeafBasisNode
387{
388 // Stores LocalFiniteElement implementations with run-time order as a function of GeometryType
389 template<typename Domain, typename Range, int dim>
390 class LagrangeRunTimeLFECache
391 {
392 public:
393 using FiniteElementType = LagrangeLocalFiniteElement<EquidistantPointSet,dim,Domain,Range>;
394
395 const FiniteElementType& get(GeometryType type)
396 {
397 auto i = data_.find(type);
398 if (i==data_.end())
399 i = data_.emplace(type,FiniteElementType(type,order_)).first;
400 return (*i).second;
401 }
402
403 std::map<GeometryType, FiniteElementType> data_;
404 unsigned int order_;
405 };
406
407 static constexpr int dim = GV::dimension;
408 static constexpr bool useDynamicOrder = (k<0);
409
410 using FiniteElementCache = typename std::conditional<(useDynamicOrder),
411 LagrangeRunTimeLFECache<typename GV::ctype, R, dim>,
412 PQkLocalFiniteElementCache<typename GV::ctype, R, dim, k>
413 >::type;
414
415public:
416
417 using size_type = std::size_t;
418 using Element = typename GV::template Codim<0>::Entity;
419 using FiniteElement = typename FiniteElementCache::FiniteElementType;
420
422 LagrangeNode() :
423 finiteElement_(nullptr),
424 element_(nullptr)
425 {}
426
428 LagrangeNode(unsigned int order) :
429 order_(order),
430 finiteElement_(nullptr),
431 element_(nullptr)
432 {
433 // Only the cache for the run-time-order case (i.e., k<0), has the 'order_' member
434 if constexpr (useDynamicOrder)
435 cache_.order_ = order;
436 }
437
439 const Element& element() const
440 {
441 return *element_;
442 }
443
448 const FiniteElement& finiteElement() const
449 {
450 return *finiteElement_;
451 }
452
454 void bind(const Element& e)
455 {
456 element_ = &e;
457 finiteElement_ = &(cache_.get(element_->type()));
458 this->setSize(finiteElement_->size());
459 }
460
461protected:
462
463 unsigned int order() const
464 {
465 return (useDynamicOrder) ? order_ : k;
466 }
467
468 // Run-time order, only valid if k<0
469 unsigned int order_;
470
471 FiniteElementCache cache_;
472 const FiniteElement* finiteElement_;
473 const Element* element_;
474};
475
476
477
478namespace BasisFactory {
479
488template<std::size_t k, typename R=double>
490{
491 return [](const auto& gridView) {
492 return LagrangePreBasis<std::decay_t<decltype(gridView)>, k, R>(gridView);
493 };
494}
495
503template<typename R=double>
504auto lagrange(int order)
505{
506 return [=](const auto& gridView) {
507 return LagrangePreBasis<std::decay_t<decltype(gridView)>, -1, R>(gridView, order);
508 };
509}
510
511} // end namespace BasisFactory
512
513
514
538template<typename GV, int k=-1, typename R=double>
540
541
542
543
544
545} // end namespace Functions
546} // end namespace Dune
547
548
549#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_LAGRANGEBASIS_HH
Global basis for given pre-basis.
Definition: defaultglobalbasis.hh:46
A pre-basis for a PQ-lagrange bases with given order.
Definition: lagrangebasis.hh:55
size_type computeDofsPerCube(std::size_t cubeDim) const
Number of degrees of freedom assigned to a cube (without the ones assigned to its faces!...
Definition: lagrangebasis.hh:350
size_type computeDofsPerSimplex(std::size_t simplexDim) const
Number of degrees of freedom assigned to a simplex (without the ones assigned to its faces!...
Definition: lagrangebasis.hh:344
LagrangeNode< GV, k, R > Node
Template mapping root tree path to type of created tree node.
Definition: lagrangebasis.hh:68
size_type dofsPerSimplex(std::size_t simplexDim) const
Number of degrees of freedom assigned to a simplex (without the ones assigned to its faces!...
Definition: lagrangebasis.hh:322
std::size_t size_type
Type used for indices and size information.
Definition: lagrangebasis.hh:65
void initializeIndices()
Initialize the global indices.
Definition: lagrangebasis.hh:96
LagrangePreBasis(const GridView &gv, unsigned int order)
Constructor for a given grid view object and run-time order.
Definition: lagrangebasis.hh:80
size_type dimension() const
Get the total dimension of the space spanned by this basis.
Definition: lagrangebasis.hh:178
LagrangePreBasis(const GridView &gv)
Constructor for a given grid view object with compile-time order.
Definition: lagrangebasis.hh:75
void update(const GridView &gv)
Update the stored grid view, to be called if the grid has changed.
Definition: lagrangebasis.hh:126
size_type size(const SizePrefix &prefix) const
Return number of possible values for next position in multi index.
Definition: lagrangebasis.hh:171
GV GridView
The grid view that the FE basis is defined on.
Definition: lagrangebasis.hh:62
const GridView & gridView() const
Obtain the grid view that the basis is defined on.
Definition: lagrangebasis.hh:120
Node makeNode() const
Create tree node.
Definition: lagrangebasis.hh:134
unsigned int order() const
Polynomial order used in the local Lagrange finite-elements.
Definition: lagrangebasis.hh:310
size_type maxNodeSize() const
Get the maximal number of DOFs associated to node for any element.
Definition: lagrangebasis.hh:184
size_type size() const
Same as size(prefix) with empty prefix.
Definition: lagrangebasis.hh:140
size_type dofsPerCube(std::size_t cubeDim) const
Number of degrees of freedom assigned to a cube (without the ones assigned to its faces!...
Definition: lagrangebasis.hh:328
Grid view abstract base class.
Definition: gridview.hh:66
Describe position of one degree of freedom.
Definition: localkey.hh:23
unsigned int index() const
Return offset within subentity.
Definition: localkey.hh:68
unsigned int codim() const
Return codim of associated entity.
Definition: localkey.hh:62
unsigned int subEntity() const
Return number of associated subentity.
Definition: localkey.hh:56
Default exception for dummy implementations.
Definition: exceptions.hh:263
Default exception class for range errors.
Definition: exceptions.hh:254
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
auto lagrange(int order)
Create a pre-basis factory that can create a Lagrange pre-basis with a run-time order.
Definition: lagrangebasis.hh:504
constexpr GeometryType prism
GeometryType representing a 3D prism.
Definition: type.hh:542
constexpr GeometryType triangle
GeometryType representing a triangle.
Definition: type.hh:518
constexpr GeometryType quadrilateral
GeometryType representing a quadrilateral (a square).
Definition: type.hh:524
constexpr GeometryType hexahedron
GeometryType representing a hexahedron.
Definition: type.hh:548
constexpr GeometryType pyramid
GeometryType representing a 3D pyramid.
Definition: type.hh:536
constexpr GeometryType tetrahedron
GeometryType representing a tetrahedron.
Definition: type.hh:530
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Convenience header that includes all implementations of Lagrange finite elements.
Dune namespace.
Definition: alignedallocator.hh:13
static constexpr T binomial(const T &n, const T &k) noexcept
calculate the binomial coefficient n over k as a constexpr
Definition: math.hh:131
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
STL namespace.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)