DUNE-FUNCTIONS (unstable)

morleybasis.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
5// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
6#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_MORLEYBASIS_HH
7#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_MORLEYBASIS_HH
8
9#include <algorithm>
10#include <type_traits>
11#include <vector>
12#include <array>
13#include <bitset>
14
15#include <dune/common/boundschecking.hh>
16#include <dune/common/exceptions.hh>
17#include <dune/common/fvector.hh>
18
19#include <dune/grid/common/mcmgmapper.hh>
20#include <dune/grid/common/rangegenerators.hh>
21
22#include <dune/localfunctions/common/localbasis.hh>
23#include <dune/localfunctions/common/localfiniteelementtraits.hh>
24#include <dune/localfunctions/common/localkey.hh>
25
26#include <dune/functions/analyticfunctions/monomialset.hh>
27#include <dune/functions/common/densevectorview.hh>
28#include <dune/functions/common/squeezetensor.hh>
30#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
31#include <dune/functions/functionspacebases/functionaldescriptor.hh>
32#include <dune/functions/functionspacebases/leafprebasismappermixin.hh>
33#include <dune/functions/functionspacebases/nodes.hh>
34#include <dune/functions/functionspacebases/transformedfiniteelementmixin.hh>
35
52namespace Dune::Functions
53{
54
55 template<class GV, class R>
56 class MorleyPreBasis;
57
75 template<class GV, class R = double>
77
78
79 namespace Impl
80 {
81
88 template<class D, class R>
89 class MorleyReferenceLocalBasis
90 {
91 public:
92 static constexpr int dim = 2;
93 using Traits = H2LocalBasisTraits<D, dim, FieldVector<D, dim>, R, 1, FieldVector<R, 1>,
94 FieldMatrix<R, 1, dim>, FieldMatrix<R, dim, dim>>;
95
96 private:
97
110 static constexpr auto getMorleyCoefficients()
111 {
112 // Define std::sqrt(2.) manually in double precision,
113 // because std::sqrt is not constexpr before C++26.
114 D sqrt2 = 0.5 * 1.414213562373095;
115
116 return Dune::FieldMatrix<D, 6,6>({{1, -1, -1, 0, 2, 0},
117 {0, 0.5, 0.5, 0.5, -1, -0.5},
118 {0, 0.5, 0.5, -0.5, -1, 0.5},
119 {0, 0, 1, 0, 0, -1},
120 {0, -1., 0, 1., 0, 0},
121 {0, sqrt2, sqrt2, -sqrt2, -2. * sqrt2, -sqrt2}});
122 }
123
124 static constexpr auto referenceBasisCoefficients = getMorleyCoefficients();
126
127 public:
128
131 static constexpr unsigned int size()
132 {
133 return 6;
134 }
135
138 static constexpr unsigned int order()
139 {
140 return 2;
141 }
142
148 void evaluateFunction(const typename Traits::DomainType &in,
149 std::vector<typename Traits::RangeType> &out) const
150 {
151 out.resize(size());
152 auto monomialValues = monomials(in);
153 multiplyWithCoefficentMatrix(referenceBasisCoefficients, monomialValues, out);
154 }
155
161 void evaluateJacobian(const typename Traits::DomainType &in,
162 std::vector<typename Traits::JacobianType> &out) const
163 {
164 out.resize(size());
165 auto monomialValues = derivative(monomials)(in);
166 multiplyWithCoefficentMatrix(referenceBasisCoefficients, monomialValues, out);
167 }
168
174 void evaluateHessian(const typename Traits::DomainType &in,
175 std::vector<typename Traits::HessianType> &out) const
176 {
177 out.resize(size());
178 auto monomialValues = derivative(derivative(monomials))(in);
179 multiplyWithCoefficentMatrix(referenceBasisCoefficients, monomialValues, out);
180 }
181
188 void partial(std::array<unsigned int, dim> order, const typename Traits::DomainType &in,
189 std::vector<typename Traits::RangeType> &out) const
190 {
191 out.resize(size());
192 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
193 if (totalOrder == 0)
194 evaluateFunction(in, out);
195 else if (totalOrder == 1)
196 {
197 evaluateJacobian(in,jacobiansBuffer_);
198 std::size_t which = std::max_element(order.begin(), order.end()) - order.begin();
199 for (auto i : Dune::range(size()))
200 out[i] = jacobiansBuffer_[i][0][which];
201 }
202 else if (totalOrder == 2)
203 {
204 evaluateHessian(in, hessianBuffer_);
205 std::size_t first, second;
206 first = std::max_element(order.begin(), order.end()) - order.begin();
207 if (order[first] == 2)
208 second = first;
209 else
210 {
211 order[first] = 0;
212 second = std::max_element(order.begin(), order.end()) - order.begin();
213 }
214 for (auto i : Dune::range(size()))
215 out[i] = hessianBuffer_[i][first][second];
216 }
217 else
218 DUNE_THROW(RangeError, "partial() not implemented for given order");
219 }
220
221 private:
222 mutable std::vector<typename Traits::JacobianType> jacobiansBuffer_;
223 mutable std::vector<typename Traits::HessianType> hessianBuffer_;
224 };
225
226
230 class MorleyLocalCoefficients
231 {
232 public:
233 using size_type = std::size_t;
234
235 MorleyLocalCoefficients()
236 {
237 for (size_type i = 0; i < 3; ++i)
238 {
239 localKeys_[i] = LocalKey(i, 2, 0); // vertex dofs
240 localKeys_[3 + i] = LocalKey(i, 1, 0); // edge dofs
241 }
242 }
243
246 static constexpr size_type size()
247 {
248 return 6;
249 }
250
253 constexpr LocalKey const &localKey(size_type i) const
254 {
255 return localKeys_[i];
256 }
257
258 private:
259 std::array<LocalKey, 6> localKeys_;
260 };
261
262
268 template<class D>
269 class MorleyLocalInterpolation
270 {
271 using size_type = std::size_t;
272 using FunctionalDescriptor = Dune::Functions::Impl::FunctionalDescriptor<2>;
273
274 public:
275
276 MorleyLocalInterpolation()
277 {
278 descriptors_[0] = FunctionalDescriptor();
279 descriptors_[1] = FunctionalDescriptor();
280 descriptors_[2] = FunctionalDescriptor();
281 descriptors_[3] = FunctionalDescriptor(1);
282 descriptors_[4] = FunctionalDescriptor(1);
283 descriptors_[5] = FunctionalDescriptor(1);
284 }
285
288 template <class Element>
289 void bind(const Element &element, const std::bitset<3> &edgeOrientation)
290 {
291 auto geometry = element.geometry();
292
293 // get global Normals and midpoints
294 auto refElement = Dune::referenceElement<double, 2>(geometry.type());
295 for (std::size_t i = 0; i < 3; ++i)
296 {
297 localVertices_[i] = refElement.position(i, 2);
298
299 localMidpoints_[i] = refElement.position(i, 1);
300 std::size_t lower = (i == 2) ? 1 : 0;
301 std::size_t upper = (i == 0) ? 1 : 2;
302
303 auto edge = geometry.global(refElement.position(upper, 2))
304 - geometry.global(refElement.position(lower, 2));
305 // normalize and orient
306 edge /= edge.two_norm() * (edgeOrientation[i] ? -1. : 1.);
307 // Rotation by pi/2. Note that Kirby rotates by 3*pi/2
308 globalNormals_[i] = {-edge[1], edge[0]};
309 }
310 }
311
319 template <class F, class C>
320 void interpolate(const F &f, std::vector<C> &out) const
321 {
322 out.resize(6);
323 auto&& df = derivative(f);
324 for (size_type i = 0; i < 3; ++i)
325 {
326 out[i] = f(localVertices_[i]);
327 out[3 + i] = squeezeTensor(df(localMidpoints_[i])).dot(globalNormals_[i]);
328 }
329 }
330
331 const FunctionalDescriptor& functionalDescriptor(size_type i) const
332 {
333 return descriptors_[i];
334 }
335
336 protected:
337 std::array<Dune::FieldVector<D, 2>, 3> globalNormals_;
338 std::array<Dune::FieldVector<D, 2>, 3> localMidpoints_;
339 std::array<Dune::FieldVector<D, 2>, 3> localVertices_;
340 std::array<FunctionalDescriptor, 6> descriptors_;
341
342 };
343
344 template<class D, class R>
345 using MorleyLocalBasisTraits = typename Impl::MorleyReferenceLocalBasis<D, R>::Traits;
346
355 template<class D, class R>
356 class MorleyLocalFiniteElement
357 : public Impl::TransformedFiniteElementMixin<MorleyLocalFiniteElement<D,R>, MorleyLocalBasisTraits<D, R>>
358 {
359 using Base = Impl::TransformedFiniteElementMixin< MorleyLocalFiniteElement<D,R>, MorleyLocalBasisTraits<D, R>>;
360 friend class Impl::TransformedLocalBasis<MorleyLocalFiniteElement<D,R>, MorleyLocalBasisTraits<D, R>>;
361 static constexpr int dim = 2;
362 public:
363
366 using size_type = std::size_t;
367 using Traits = LocalFiniteElementTraits<
368 Impl::TransformedLocalBasis<MorleyLocalFiniteElement<D,R>, MorleyLocalBasisTraits<D, R>>,
369 Impl::MorleyLocalCoefficients,
370 Impl::MorleyLocalInterpolation<D>>;
371
375 const typename Traits::LocalCoefficientsType &localCoefficients() const
376 {
377 return coefficients_;
378 }
379
382 const typename Traits::LocalInterpolationType &localInterpolation() const
383 {
384 return interpolation_;
385 }
386
389 static constexpr GeometryType type()
390 {
391 return GeometryTypes::simplex(dim);
392 }
393
396 static constexpr size_type size()
397 {
398 return 6;
399 }
400
403 template<class Element>
404 void bind(std::bitset<3> const& data, Element const &e)
405 {
406 edgeOrientation_ = data;
407
408 fillMatrix(e.geometry());
409 interpolation_.bind(e, edgeOrientation_);
410 }
411
412 protected:
413
416 Impl::MorleyReferenceLocalBasis<D, R> const& referenceLocalBasis() const
417 {
418 return basis_;
419 }
420
426 template<class InputValues, class OutputValues>
427 void transform(InputValues const& inValues, OutputValues& outValues) const
428 {
429 // Here we cannot directly use
430 // mat_.mv(inValues, outValues);
431 // because mv expects the DenseVector interface.
432 auto inValuesDenseVector = Impl::DenseVectorView(inValues);
433 auto outValuesDenseVector = Impl::DenseVectorView(outValues);
434 mat_.mv(inValuesDenseVector, outValuesDenseVector);
435 }
436
437 private:
438
445 template<class Geometry>
446 void fillMatrix(Geometry const &geometry)
447 {
448 std::array<R, 3> B_11;
449 std::array<R, 3> B_12;
450 std::array<R, 3> l_inv;
451
452 std::array<Dune::FieldVector<R, 2>, 3> referenceTangents;
453 std::array<Dune::FieldVector<R, 2>, 3> globalTangents;
454
455 // By default, edges point from the vertex with the smaller index
456 // to the vertex with the larger index.
457
458 // get local and global Tangents
459 auto refElement = Dune::referenceElement<double, 2>(geometry.type());
460 auto x = refElement.position(0,0);
461 for (std::size_t i = 0; i < 3; ++i)
462 {
463 std::size_t lower = (i == 2) ? 1 : 0;
464 std::size_t upper = (i == 0) ? 1 : 2;
465 auto edge = refElement.position(upper, 2) - refElement.position(lower, 2);
466
467 referenceTangents[i] = edge / edge.two_norm();
468
469 auto globalEdge = geometry.global(refElement.position(upper, 2))
470 - geometry.global(refElement.position(lower, 2));
471
472 l_inv[i] = 1. / globalEdge.two_norm();
473 globalTangents[i] = globalEdge * l_inv[i];
474 }
475
476 auto jacobianTransposed = geometry.jacobianTransposed(x);
477 for (std::size_t i = 0; i < 3; ++i)
478 {
479 B_11[i] = -referenceTangents[i][1]
480 *(-globalTangents[i][1] * jacobianTransposed[0][0]
481 + globalTangents[i][0] * jacobianTransposed[0][1])
482 + referenceTangents[i][0]
483 *(-globalTangents[i][1] * jacobianTransposed[1][0]
484 + globalTangents[i][0] * jacobianTransposed[1][1]);
485 B_12[i] = -referenceTangents[i][1]
486 *(globalTangents[i][0] * jacobianTransposed[0][0]
487 + globalTangents[i][1] * jacobianTransposed[0][1])
488 + referenceTangents[i][0]
489 *(globalTangents[i][0] * jacobianTransposed[1][0]
490 + globalTangents[i][1] * jacobianTransposed[1][1]);
491 }
492
493 // Actually setup matrix
494 int sign = -1;
495 for (std::size_t i = 0; i < 3; ++i)
496 {
497 mat_[i][i] = 1.;
498 for (std::size_t j = 0; j < 3; ++j)
499 {
500 if (j != (2 - i)) // dune specific edge order
501 {
502 mat_[j][3 + i] = sign * B_12[i] * l_inv[i];
503 sign *= -1;
504 }
505 }
506 mat_[3 + i][3 + i] = (edgeOrientation_[i] ? -1. : 1.) * B_11[i];
507 }
508 }
509
510 // a finite element consists of a basis, coeffiecents and an interpolation
511 typename Impl::MorleyReferenceLocalBasis<D, R> basis_;
512 typename Traits::LocalCoefficientsType coefficients_;
513 typename Traits::LocalInterpolationType interpolation_;
514 // This is the matrix M in Kirbys paper
515 Dune::FieldMatrix<R, 6, 6>mat_;
516 // the local state, i.e. a collection of global information restricted to this element
517 std::bitset<3> edgeOrientation_;
518
519 };
520
521 } // end namespace Impl
522
523
524
525 // *****************************************************************************
526 // This is the reusable part of the basis. It contains
527 //
528 // MorleyPreBasis
529 // MorleyNode
530 //
531 // The pre-basis allows to create the others and is the owner of possible shared
532 // state. These components do _not_ depend on the global basis and local view
533 // and can be used without a global basis.
534 // *****************************************************************************
535
536 template<class GV, class R>
537 class MorleyNode
538 : public LeafBasisNode
539 {
540 using Mapper = Dune::MultipleCodimMultipleGeomTypeMapper<GV>;
541 public:
542 using size_type = std::size_t;
543 using Element = typename GV::template Codim<0>::Entity;
544 using FiniteElement = typename Impl::MorleyLocalFiniteElement<typename GV::ctype, R>;
545
546
547 MorleyNode(Mapper const& m, std::vector<std::bitset<3>> const& data)
548 : mapper_(&m), data_(&data)
549 {}
550
552 Element const &element() const
553 {
554 return *element_;
555 }
556
562 FiniteElement const &finiteElement() const
563 {
564 return finiteElement_;
565 }
566
568 void bind(Element const &e)
569 {
570 element_ = &e;
571 finiteElement_.bind((*data_)[mapper_->index(e)], *element_);
572 this->setSize(finiteElement_.size());
573 }
574
576 unsigned int order() const
577 {
578 return finiteElement_.localBasis().order();
579 }
580
581 protected:
582 FiniteElement finiteElement_;
583 Element const* element_;
584 Mapper const* mapper_;
585 std::vector<std::bitset<3>> const* data_;
586 };
587
588
598 template<class GV, class R>
600 : public LeafPreBasisMapperMixin<GV>
601 {
603 using SubEntityMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GV>;
604 using Element = typename GV::template Codim<0>::Entity;
605 using D = typename GV::ctype;
606 static const std::size_t dim = GV::dimension;
607
608 // helper methods to assign each subentity the number of dofs. Used by the LeafPreBasisMapperMixin.
609 static constexpr auto morleyMapperLayout(Dune::GeometryType type, int gridDim)
610 {
611 assert(gridDim == 2);
612 if (type.isVertex())
613 return 1; // one evaluation dof per vertex
614 if (type.isLine())
615 return 1;
616 if ((type.isTriangle()) )
617 return 0;
618 else
619 return 0;
620 }
621
622 public:
624 using GridView = GV;
625
627 using size_type = std::size_t;
628
630 using Node = MorleyNode<GridView, R>;
631
633 MorleyPreBasis(const GV &gv)
634 : Base(gv, morleyMapperLayout)
635 , mapper_({gv, mcmgElementLayout()})
636 {
637 data_ = Impl::computeEdgeOrientations(mapper_);
638 }
639
641 void update(GridView const &gv)
642 {
643 Base::update(gv);
644 data_ = Impl::computeEdgeOrientations(mapper_);
645 }
646
651 {
652 return Node{mapper_, data_};
653 }
654
655 protected:
656
657 SubEntityMapper mapper_;
658 std::vector<std::bitset<3>> data_;
659
660 }; // class MorleyPreBasis
661
662 namespace BasisFactory
663 {
664
671 template<class R = double>
672 auto morley()
673 {
674 return [=](auto const &gridView) {
675 return MorleyPreBasis<std::decay_t<decltype(gridView)>, R>(gridView);
676 };
677 }
678
679 } // end namespace BasisFactory
680
681} // end namespace Dune::Functions
682
683#endif
Global basis for given pre-basis.
Definition: defaultglobalbasis.hh:53
A generic MixIn class for PreBasis with flat indices computed from a mapper.
Definition: leafprebasismappermixin.hh:62
std::size_t size_type
Type used for index digits.
Definition: leafprebasismappermixin.hh:71
void update(const GridView &gv)
Update the stored GridView.
Definition: leafprebasismappermixin.hh:101
GV GridView
Type of the associated GridView.
Definition: leafprebasismappermixin.hh:68
A pre-basis for a Morleybasis.
Definition: morleybasis.hh:601
void update(GridView const &gv)
Update the stored grid view, to be called if the grid has changed.
Definition: morleybasis.hh:641
Node makeNode() const
Create tree node.
Definition: morleybasis.hh:650
MorleyNode< GridView, R > Node
Template mapping root tree path to type of created tree node.
Definition: morleybasis.hh:630
MorleyPreBasis(const GV &gv)
Constructor for a given grid view object.
Definition: morleybasis.hh:633
This file provides an implementation of the cubic Hermite finite element in 1 to 3 dimensions.
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:43
auto morley()
construct a PreBasisFactory for the Morley Finite Element
Definition: morleybasis.hh:672
Function, which evaluates all monomials up to degree maxDegree in a given coordinate.
Definition: monomialset.hh:64
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Jan 9, 23:34, 2026)