Dune Core Modules (2.9.0)

lagrangecube.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
7
8#include <array>
9#include <numeric>
10
13#include <dune/common/math.hh>
14
15#include <dune/geometry/referenceelements.hh>
16
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19#include <dune/localfunctions/common/localinterpolation.hh>
20#include <dune/localfunctions/common/localkey.hh>
21
22namespace Dune { namespace Impl
23{
24 // Forward declaration
25 template<class LocalBasis>
26 class LagrangeCubeLocalInterpolation;
27
38 template<class D, class R, unsigned int dim, unsigned int k>
39 class LagrangeCubeLocalBasis
40 {
41 friend class LagrangeCubeLocalInterpolation<LagrangeCubeLocalBasis<D,R,dim,k> >;
42
43 // i-th Lagrange polynomial of degree k in one dimension
44 static R p(unsigned int i, D x)
45 {
46 R result(1.0);
47 for (unsigned int j=0; j<=k; j++)
48 if (j!=i) result *= (k*x-j)/((int)i-(int)j);
49 return result;
50 }
51
52 // derivative of ith Lagrange polynomial of degree k in one dimension
53 static R dp(unsigned int i, D x)
54 {
55 R result(0.0);
56
57 for (unsigned int j=0; j<=k; j++)
58 {
59 if (j!=i)
60 {
61 R prod( (k*1.0)/((int)i-(int)j) );
62 for (unsigned int l=0; l<=k; l++)
63 if (l!=i && l!=j)
64 prod *= (k*x-l)/((int)i-(int)l);
65 result += prod;
66 }
67 }
68 return result;
69 }
70
71 // Second derivative of j-th Lagrange polynomial of degree k in one dimension
72 // Formula and notation taken from https://en.wikipedia.org/wiki/Lagrange_polynomial#Derivatives
73 static R ddp(unsigned int j, D x)
74 {
75 R result(0.0);
76
77 for (unsigned int i=0; i<=k; i++)
78 {
79 if (i==j)
80 continue;
81
82 R sum(0);
83
84 for (unsigned int m=0; m<=k; m++)
85 {
86 if (m==i || m==j)
87 continue;
88
89 R prod( (k*1.0)/((int)j-(int)m) );
90 for (unsigned int l=0; l<=k; l++)
91 if (l!=i && l!=j && l!=m)
92 prod *= (k*x-l)/((int)j-(int)l);
93 sum += prod;
94 }
95
96 result += sum * ( (k*1.0)/((int)j-(int)i) );
97 }
98
99 return result;
100 }
101
102 // Return i as a d-digit number in the (k+1)-nary system
103 static std::array<unsigned int,dim> multiindex (unsigned int i)
104 {
105 std::array<unsigned int,dim> alpha;
106 for (unsigned int j=0; j<dim; j++)
107 {
108 alpha[j] = i % (k+1);
109 i = i/(k+1);
110 }
111 return alpha;
112 }
113
114 public:
115 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
116
119 static constexpr unsigned int size ()
120 {
121 return power(k+1, dim);
122 }
123
125 void evaluateFunction(const typename Traits::DomainType& x,
126 std::vector<typename Traits::RangeType>& out) const
127 {
128 out.resize(size());
129
130 // Specialization for zero-order case
131 if (k==0)
132 {
133 out[0] = 1;
134 return;
135 }
136
137 if (k==1)
138 {
139 for (size_t i=0; i<size(); i++)
140 {
141 out[i] = 1;
142
143 for (unsigned int j=0; j<dim; j++)
144 // if j-th bit of i is set multiply with x[j], else with 1-x[j]
145 out[i] *= (i & (1<<j)) ? x[j] : 1-x[j];
146 }
147 return;
148 }
149
150 // General case
151 for (size_t i=0; i<size(); i++)
152 {
153 // convert index i to multiindex
154 std::array<unsigned int,dim> alpha(multiindex(i));
155
156 // initialize product
157 out[i] = 1.0;
158
159 // dimension by dimension
160 for (unsigned int j=0; j<dim; j++)
161 out[i] *= p(alpha[j],x[j]);
162 }
163 }
164
170 void evaluateJacobian(const typename Traits::DomainType& x,
171 std::vector<typename Traits::JacobianType>& out) const
172 {
173 out.resize(size());
174
175 // Specialization for k==0
176 if (k==0)
177 {
178 std::fill(out[0][0].begin(), out[0][0].end(), 0);
179 return;
180 }
181
182 // Specialization for k==1
183 if (k==1)
184 {
185 // Loop over all shape functions
186 for (size_t i=0; i<size(); i++)
187 {
188 // Loop over all coordinate directions
189 for (unsigned int j=0; j<dim; j++)
190 {
191 // Initialize: the overall expression is a product
192 // if j-th bit of i is set to 1, else -11
193 out[i][0][j] = (i & (1<<j)) ? 1 : -1;
194
195 for (unsigned int l=0; l<dim; l++)
196 {
197 if (j!=l)
198 // if l-th bit of i is set multiply with x[l], else with 1-x[l]
199 out[i][0][j] *= (i & (1<<l)) ? x[l] : 1-x[l];
200 }
201 }
202 }
203 return;
204 }
205
206 // The general case
207
208 // Loop over all shape functions
209 for (size_t i=0; i<size(); i++)
210 {
211 // convert index i to multiindex
212 std::array<unsigned int,dim> alpha(multiindex(i));
213
214 // Loop over all coordinate directions
215 for (unsigned int j=0; j<dim; j++)
216 {
217 // Initialize: the overall expression is a product
218 // if j-th bit of i is set to -1, else 1
219 out[i][0][j] = dp(alpha[j],x[j]);
220
221 // rest of the product
222 for (unsigned int l=0; l<dim; l++)
223 if (l!=j)
224 out[i][0][j] *= p(alpha[l],x[l]);
225 }
226 }
227 }
228
235 void partial(const std::array<unsigned int,dim>& order,
236 const typename Traits::DomainType& in,
237 std::vector<typename Traits::RangeType>& out) const
238 {
239 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
240
241 out.resize(size());
242
243 if (k==0)
244 {
245 out[0] = (totalOrder==0);
246 return;
247 }
248
249 if (k==1)
250 {
251 if (totalOrder == 0)
252 {
253 evaluateFunction(in, out);
254 }
255 else if (totalOrder == 1)
256 {
257 out.resize(size());
258
259 auto direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
260 if (direction >= dim)
261 DUNE_THROW(RangeError, "Direction of partial derivative not found!");
262
263 // Loop over all shape functions
264 for (std::size_t i = 0; i < size(); ++i)
265 {
266 // Initialize: the overall expression is a product
267 // if j-th bit of i is set to 1, otherwise to -1
268 out[i] = (i & (1<<direction)) ? 1 : -1;
269
270 for (unsigned int j = 0; j < dim; ++j)
271 {
272 if (direction != j)
273 // if j-th bit of i is set multiply with in[j], else with 1-in[j]
274 out[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
275 }
276 }
277 }
278 else if (totalOrder == 2)
279 {
280
281 for (size_t i=0; i<size(); i++)
282 {
283 // convert index i to multiindex
284 std::array<unsigned int,dim> alpha(multiindex(i));
285
286 // Initialize: the overall expression is a product
287 out[i][0] = 1.0;
288
289 // rest of the product
290 for (std::size_t l=0; l<dim; l++)
291 {
292 switch (order[l])
293 {
294 case 0:
295 out[i][0] *= p(alpha[l],in[l]);
296 break;
297 case 1:
298 //std::cout << "dp: " << dp(alpha[l],in[l]) << std::endl;
299 out[i][0] *= dp(alpha[l],in[l]);
300 break;
301 case 2:
302 out[i][0] *= 0;
303 break;
304 default:
305 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
306 }
307 }
308 }
309 }
310 else
311 DUNE_THROW(NotImplemented, "Partial derivative of order " << totalOrder << " is not implemented!");
312
313 return;
314 }
315
316 // The case k>1
317
318 // Loop over all shape functions
319 for (size_t i=0; i<size(); i++)
320 {
321 // convert index i to multiindex
322 std::array<unsigned int,dim> alpha(multiindex(i));
323
324 // Initialize: the overall expression is a product
325 out[i][0] = 1.0;
326
327 // rest of the product
328 for (std::size_t l=0; l<dim; l++)
329 {
330 switch (order[l])
331 {
332 case 0:
333 out[i][0] *= p(alpha[l],in[l]);
334 break;
335 case 1:
336 out[i][0] *= dp(alpha[l],in[l]);
337 break;
338 case 2:
339 out[i][0] *= ddp(alpha[l],in[l]);
340 break;
341 default:
342 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
343 }
344 }
345 }
346 }
347
349 static constexpr unsigned int order ()
350 {
351 return k;
352 }
353 };
354
360 template<unsigned int dim, unsigned int k>
361 class LagrangeCubeLocalCoefficients
362 {
363 // Return i as a d-digit number in the (k+1)-nary system
364 static std::array<unsigned int,dim> multiindex (unsigned int i)
365 {
366 std::array<unsigned int,dim> alpha;
367 for (unsigned int j=0; j<dim; j++)
368 {
369 alpha[j] = i % (k+1);
370 i = i/(k+1);
371 }
372 return alpha;
373 }
374
376 void setup1d(std::vector<unsigned int>& subEntity)
377 {
378 assert(k>0);
379
380 unsigned lastIndex=0;
381
382 /* edge and vertex numbering
383 0----0----1
384 */
385
386 // edge (0)
387 subEntity[lastIndex++] = 0; // corner 0
388 for (unsigned i = 0; i < k - 1; ++i)
389 subEntity[lastIndex++] = 0; // inner dofs of element (0)
390
391 subEntity[lastIndex++] = 1; // corner 1
392
393 assert(power(k+1,dim)==lastIndex);
394 }
395
396 void setup2d(std::vector<unsigned int>& subEntity)
397 {
398 assert(k>0);
399
400 unsigned lastIndex=0;
401
402 // LocalKey: entity number, entity codim, dof indices within each entity
403 /* edge and vertex numbering
404 2----3----3
405 | |
406 | |
407 0 1
408 | |
409 | |
410 0----2----1
411 */
412
413 // lower edge (2)
414 subEntity[lastIndex++] = 0; // corner 0
415 for (unsigned i = 0; i < k - 1; ++i)
416 subEntity[lastIndex++] = 2; // inner dofs of lower edge (2)
417
418 subEntity[lastIndex++] = 1; // corner 1
419
420 // iterate from bottom to top over inner edge dofs
421 for (unsigned e = 0; e < k - 1; ++e) {
422 subEntity[lastIndex++] = 0; // left edge (0)
423 for (unsigned i = 0; i < k - 1; ++i)
424 subEntity[lastIndex++] = 0; // face dofs
425 subEntity[lastIndex++] = 1; // right edge (1)
426 }
427
428 // upper edge (3)
429 subEntity[lastIndex++] = 2; // corner 2
430 for (unsigned i = 0; i < k - 1; ++i)
431 subEntity[lastIndex++] = 3; // inner dofs of upper edge (3)
432
433 subEntity[lastIndex++] = 3; // corner 3
434
435 assert(power(k+1,dim)==lastIndex);
436 }
437
438 void setup3d(std::vector<unsigned int>& subEntity)
439 {
440 assert(k>0);
441
442 unsigned lastIndex=0;
443#ifndef NDEBUG
444 const unsigned numIndices = power(k+1,dim);
445 const unsigned numFaceIndices = power(k+1,dim-1);
446#endif
447 const unsigned numInnerEdgeDofs = k-1;
448
449 // LocalKey: entity number, entity codim, dof indices within each entity
450 /* edge and vertex numbering
451
452 6---(11)--7 6---------7
453 /| /| /| (5) /|
454 (8)| (9)| / | top / |
455 / (2) / (3) / |(3)bac/k |
456 4---(10)--5 | 4---------5 |
457 | | | | left|(0)| |(1)|right
458 | 2--(7)|---3 | 2-----|---3
459 (0) / (1) / |(2)front | /
460 |(4) |(5) | / (4) | /
461 |/ |/ |/ bottom |/
462 0---(6)---1 0---------1
463 */
464
465 // bottom face (4)
466 lastIndex=0;
467 // lower edge (6)
468 subEntity[lastIndex++] = 0; // corner 0
469 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
470 subEntity[lastIndex++] = 6; // inner dofs of lower edge (6)
471
472 subEntity[lastIndex++] = 1; // corner 1
473
474 // iterate from bottom to top over inner edge dofs
475 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
476 subEntity[lastIndex++] = 4; // left edge (4)
477 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
478 subEntity[lastIndex++] = 4; // inner face dofs
479 subEntity[lastIndex++] = 5; // right edge (5)
480 }
481
482 // upper edge (7)
483 subEntity[lastIndex++] = 2; // corner 2
484 for (unsigned i = 0; i < k - 1; ++i)
485 subEntity[lastIndex++] = 7; // inner dofs of upper edge (7)
486 subEntity[lastIndex++] = 3; // corner 3
487
488 assert(numFaceIndices==lastIndex); // added 1 face so far
490
492 for(unsigned f = 0; f < numInnerEdgeDofs; ++f) {
493
494 // lower edge (connecting edges 0 and 1)
495 subEntity[lastIndex++] = 0; // dof on edge 0
496 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
497 subEntity[lastIndex++] = 2; // dof in front face
498 subEntity[lastIndex++] = 1; // dof on edge 1
499
500 // iterate from bottom to top over inner edge dofs
501 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
502 subEntity[lastIndex++] = 0; // on left face (0)
503 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
504 subEntity[lastIndex++] = 0; // volume dofs
505 subEntity[lastIndex++] = 1; // right face (1)
506 }
507
508 // upper edge (connecting edges 0 and 1)
509 subEntity[lastIndex++] = 2; // dof on edge 2
510 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
511 subEntity[lastIndex++] = 3; // dof on rear face (3)
512 subEntity[lastIndex++] = 3; // dof on edge 3
513
514 assert(lastIndex==(f+1+1)*numFaceIndices);
515 }
516
518 // lower edge (10)
519 subEntity[lastIndex++] = 4; // corner 4
520 for (unsigned i = 0; i < k - 1; ++i)
521 subEntity[lastIndex++] = 10; // inner dofs on lower edge (10)
522 subEntity[lastIndex++] = 5; // corner 5
523
524 // iterate from bottom to top over inner edge dofs
525 for (unsigned e = 0; e < k - 1; ++e) {
526 subEntity[lastIndex++] = 8; // left edge (8)
527 for (unsigned i = 0; i < k - 1; ++i)
528 subEntity[lastIndex++] = 5; // face dofs
529 subEntity[lastIndex++] = 9; // right edge (9)
530 }
531
532 // upper edge (11)
533 subEntity[lastIndex++] = 6; // corner 6
534 for (unsigned i = 0; i < k - 1; ++i)
535 subEntity[lastIndex++] = 11; // inner dofs of upper edge (11)
536 subEntity[lastIndex++] = 7; // corner 7
537
538 assert(numIndices==lastIndex);
539 }
540
541 public:
543 LagrangeCubeLocalCoefficients ()
544 : localKeys_(size())
545 {
546 if (k==0)
547 {
548 localKeys_[0] = LocalKey(0,0,0);
549 return;
550 }
551
552 if (k==1)
553 {
554 for (std::size_t i=0; i<size(); i++)
555 localKeys_[i] = LocalKey(i,dim,0);
556 return;
557 }
558
559 // Now: the general case
560
561 // Set up array of codimension-per-dof-number
562 std::vector<unsigned int> codim(size());
563
564 for (std::size_t i=0; i<codim.size(); i++)
565 {
566 codim[i] = 0;
567
568 // Codimension gets increased by 1 for each coordinate direction
569 // where dof is on boundary
570 std::array<unsigned int,dim> mIdx = multiindex(i);
571 for (unsigned int j=0; j<dim; j++)
572 if (mIdx[j]==0 or mIdx[j]==k)
573 codim[i]++;
574 }
575
576 // Set up index vector (the index of the dof in the set of dofs of a given subentity)
577 // Algorithm: the 'index' has the same ordering as the dof number 'i'.
578 // To make it consecutive we interpret 'i' in the (k+1)-adic system, omit all digits
579 // that correspond to axes where the dof is on the element boundary, and transform the
580 // rest to the (k-1)-adic system.
581 std::vector<unsigned int> index(size());
582
583 for (std::size_t i=0; i<size(); i++)
584 {
585 index[i] = 0;
586
587 std::array<unsigned int,dim> mIdx = multiindex(i);
588
589 for (int j=dim-1; j>=0; j--)
590 if (mIdx[j]>0 && mIdx[j]<k)
591 index[i] = (k-1)*index[i] + (mIdx[j]-1);
592 }
593
594 // Set up entity and dof numbers for each (supported) dimension separately
595 std::vector<unsigned int> subEntity(size());
596
597 if (dim==1) {
598
599 setup1d(subEntity);
600
601 } else if (dim==2) {
602
603 setup2d(subEntity);
604
605 } else if (dim==3) {
606
607 setup3d(subEntity);
608
609 } else
610 DUNE_THROW(Dune::NotImplemented, "LagrangeCubeLocalCoefficients for order " << k << " and dim == " << dim);
611
612 for (size_t i=0; i<size(); i++)
613 localKeys_[i] = LocalKey(subEntity[i], codim[i], index[i]);
614 }
615
617 static constexpr std::size_t size ()
618 {
619 return power(k+1,dim);
620 }
621
623 const LocalKey& localKey (std::size_t i) const
624 {
625 return localKeys_[i];
626 }
627
628 private:
629 std::vector<LocalKey> localKeys_;
630 };
631
636 template<class LocalBasis>
637 class LagrangeCubeLocalInterpolation
638 {
639 public:
640
648 template<typename F, typename C>
649 void interpolate (const F& ff, std::vector<C>& out) const
650 {
651 constexpr auto dim = LocalBasis::Traits::dimDomain;
652 constexpr auto k = LocalBasis::order();
653 using D = typename LocalBasis::Traits::DomainFieldType;
654
655 typename LocalBasis::Traits::DomainType x;
656 auto&& f = Impl::makeFunctionWithCallOperator<typename LocalBasis::Traits::DomainType>(ff);
657
658 out.resize(LocalBasis::size());
659
660 // Specialization for zero-order case
661 if (k==0)
662 {
663 auto center = ReferenceElements<D,dim>::cube().position(0,0);
664 out[0] = f(center);
665 return;
666 }
667
668 // Specialization for first-order case
669 if (k==1)
670 {
671 for (unsigned int i=0; i<LocalBasis::size(); i++)
672 {
673 // Generate coordinate of the i-th corner of the reference cube
674 for (int j=0; j<dim; j++)
675 x[j] = (i & (1<<j)) ? 1.0 : 0.0;
676
677 out[i] = f(x);
678 }
679 return;
680 }
681
682 // The general case
683 for (unsigned int i=0; i<LocalBasis::size(); i++)
684 {
685 // convert index i to multiindex
686 std::array<unsigned int,dim> alpha(LocalBasis::multiindex(i));
687
688 // Generate coordinate of the i-th Lagrange point
689 for (unsigned int j=0; j<dim; j++)
690 x[j] = (1.0*alpha[j])/k;
691
692 out[i] = f(x);
693 }
694 }
695
696 };
697
698} } // namespace Dune::Impl
699
700namespace Dune
701{
709 template<class D, class R, int dim, int k>
711 {
712 public:
716 Impl::LagrangeCubeLocalCoefficients<dim,k>,
717 Impl::LagrangeCubeLocalInterpolation<Impl::LagrangeCubeLocalBasis<D,R,dim,k> > >;
718
725
728 const typename Traits::LocalBasisType& localBasis () const
729 {
730 return basis_;
731 }
732
736 {
737 return coefficients_;
738 }
739
743 {
744 return interpolation_;
745 }
746
748 static constexpr std::size_t size ()
749 {
750 return power(k+1,dim);
751 }
752
755 static constexpr GeometryType type ()
756 {
757 return GeometryTypes::cube(dim);
758 }
759
760 private:
761 Impl::LagrangeCubeLocalBasis<D,R,dim,k> basis_;
762 Impl::LagrangeCubeLocalCoefficients<dim,k> coefficients_;
763 Impl::LagrangeCubeLocalInterpolation<Impl::LagrangeCubeLocalBasis<D,R,dim,k> > interpolation_;
764 };
765
766} // namespace Dune
767
768#endif // DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
Lagrange finite element for cubes with arbitrary compile-time dimension and polynomial order.
Definition: lagrangecube.hh:711
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangecube.hh:728
LagrangeCubeLocalFiniteElement()
Default constructor.
Definition: lagrangecube.hh:724
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangecube.hh:742
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangecube.hh:755
static constexpr std::size_t size()
The number of shape functions.
Definition: lagrangecube.hh:748
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangecube.hh:735
Default exception for dummy implementations.
Definition: exceptions.hh:263
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:472
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
static const ReferenceElement & cube()
get hypercube reference elements
Definition: referenceelements.hh:210
D DomainType
domain type
Definition: localbasis.hh:42
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)