6#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_MORLEYBASIS_HH
7#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_MORLEYBASIS_HH
15#include <dune/common/boundschecking.hh>
16#include <dune/common/exceptions.hh>
17#include <dune/common/fvector.hh>
19#include <dune/grid/common/mcmgmapper.hh>
20#include <dune/grid/common/rangegenerators.hh>
22#include <dune/localfunctions/common/localbasis.hh>
23#include <dune/localfunctions/common/localfiniteelementtraits.hh>
24#include <dune/localfunctions/common/localkey.hh>
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>
52namespace Dune::Functions
55 template<
class GV,
class R>
75 template<
class GV,
class R =
double>
88 template<
class D,
class R>
89 class MorleyReferenceLocalBasis
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>>;
110 static constexpr auto getMorleyCoefficients()
114 D sqrt2 = 0.5 * 1.414213562373095;
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},
120 {0, -1., 0, 1., 0, 0},
121 {0, sqrt2, sqrt2, -sqrt2, -2. * sqrt2, -sqrt2}});
124 static constexpr auto referenceBasisCoefficients = getMorleyCoefficients();
131 static constexpr unsigned int size()
138 static constexpr unsigned int order()
148 void evaluateFunction(
const typename Traits::DomainType &in,
149 std::vector<typename Traits::RangeType> &out)
const
152 auto monomialValues = monomials(in);
153 multiplyWithCoefficentMatrix(referenceBasisCoefficients, monomialValues, out);
161 void evaluateJacobian(
const typename Traits::DomainType &in,
162 std::vector<typename Traits::JacobianType> &out)
const
165 auto monomialValues =
derivative(monomials)(in);
166 multiplyWithCoefficentMatrix(referenceBasisCoefficients, monomialValues, out);
174 void evaluateHessian(
const typename Traits::DomainType &in,
175 std::vector<typename Traits::HessianType> &out)
const
179 multiplyWithCoefficentMatrix(referenceBasisCoefficients, monomialValues, out);
188 void partial(std::array<unsigned int, dim> order,
const typename Traits::DomainType &in,
189 std::vector<typename Traits::RangeType> &out)
const
192 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
194 evaluateFunction(in, out);
195 else if (totalOrder == 1)
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];
202 else if (totalOrder == 2)
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)
212 second = std::max_element(order.begin(), order.end()) - order.begin();
214 for (
auto i : Dune::range(size()))
215 out[i] = hessianBuffer_[i][first][second];
218 DUNE_THROW(RangeError,
"partial() not implemented for given order");
222 mutable std::vector<typename Traits::JacobianType> jacobiansBuffer_;
223 mutable std::vector<typename Traits::HessianType> hessianBuffer_;
230 class MorleyLocalCoefficients
233 using size_type = std::size_t;
235 MorleyLocalCoefficients()
237 for (size_type i = 0; i < 3; ++i)
239 localKeys_[i] = LocalKey(i, 2, 0);
240 localKeys_[3 + i] = LocalKey(i, 1, 0);
246 static constexpr size_type size()
253 constexpr LocalKey
const &localKey(size_type i)
const
255 return localKeys_[i];
259 std::array<LocalKey, 6> localKeys_;
269 class MorleyLocalInterpolation
271 using size_type = std::size_t;
272 using FunctionalDescriptor = Dune::Functions::Impl::FunctionalDescriptor<2>;
276 MorleyLocalInterpolation()
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);
288 template <
class Element>
289 void bind(
const Element &element,
const std::bitset<3> &edgeOrientation)
291 auto geometry = element.geometry();
294 auto refElement = Dune::referenceElement<double, 2>(geometry.type());
295 for (std::size_t i = 0; i < 3; ++i)
297 localVertices_[i] = refElement.position(i, 2);
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;
303 auto edge = geometry.global(refElement.position(upper, 2))
304 - geometry.global(refElement.position(lower, 2));
306 edge /= edge.two_norm() * (edgeOrientation[i] ? -1. : 1.);
308 globalNormals_[i] = {-edge[1], edge[0]};
319 template <
class F,
class C>
320 void interpolate(
const F &f, std::vector<C> &out)
const
324 for (size_type i = 0; i < 3; ++i)
326 out[i] = f(localVertices_[i]);
327 out[3 + i] = squeezeTensor(df(localMidpoints_[i])).dot(globalNormals_[i]);
331 const FunctionalDescriptor& functionalDescriptor(size_type i)
const
333 return descriptors_[i];
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_;
344 template<
class D,
class R>
345 using MorleyLocalBasisTraits =
typename Impl::MorleyReferenceLocalBasis<D, R>::Traits;
355 template<
class D,
class R>
356 class MorleyLocalFiniteElement
357 :
public Impl::TransformedFiniteElementMixin<MorleyLocalFiniteElement<D,R>, MorleyLocalBasisTraits<D, R>>
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;
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>>;
375 const typename Traits::LocalCoefficientsType &localCoefficients()
const
377 return coefficients_;
382 const typename Traits::LocalInterpolationType &localInterpolation()
const
384 return interpolation_;
389 static constexpr GeometryType type()
391 return GeometryTypes::simplex(dim);
396 static constexpr size_type size()
403 template<
class Element>
404 void bind(std::bitset<3>
const& data, Element
const &e)
406 edgeOrientation_ = data;
408 fillMatrix(e.geometry());
409 interpolation_.bind(e, edgeOrientation_);
416 Impl::MorleyReferenceLocalBasis<D, R>
const& referenceLocalBasis()
const
426 template<
class InputValues,
class OutputValues>
427 void transform(InputValues
const& inValues, OutputValues& outValues)
const
432 auto inValuesDenseVector = Impl::DenseVectorView(inValues);
433 auto outValuesDenseVector = Impl::DenseVectorView(outValues);
434 mat_.mv(inValuesDenseVector, outValuesDenseVector);
445 template<
class Geometry>
446 void fillMatrix(Geometry
const &geometry)
448 std::array<R, 3> B_11;
449 std::array<R, 3> B_12;
450 std::array<R, 3> l_inv;
452 std::array<Dune::FieldVector<R, 2>, 3> referenceTangents;
453 std::array<Dune::FieldVector<R, 2>, 3> globalTangents;
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)
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);
467 referenceTangents[i] = edge / edge.two_norm();
469 auto globalEdge = geometry.global(refElement.position(upper, 2))
470 - geometry.global(refElement.position(lower, 2));
472 l_inv[i] = 1. / globalEdge.two_norm();
473 globalTangents[i] = globalEdge * l_inv[i];
476 auto jacobianTransposed = geometry.jacobianTransposed(x);
477 for (std::size_t i = 0; i < 3; ++i)
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]);
495 for (std::size_t i = 0; i < 3; ++i)
498 for (std::size_t j = 0; j < 3; ++j)
502 mat_[j][3 + i] = sign * B_12[i] * l_inv[i];
506 mat_[3 + i][3 + i] = (edgeOrientation_[i] ? -1. : 1.) * B_11[i];
511 typename Impl::MorleyReferenceLocalBasis<D, R> basis_;
512 typename Traits::LocalCoefficientsType coefficients_;
513 typename Traits::LocalInterpolationType interpolation_;
515 Dune::FieldMatrix<R, 6, 6>mat_;
517 std::bitset<3> edgeOrientation_;
536 template<
class GV,
class R>
538 :
public LeafBasisNode
540 using Mapper = Dune::MultipleCodimMultipleGeomTypeMapper<GV>;
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>;
547 MorleyNode(Mapper
const& m, std::vector<std::bitset<3>>
const& data)
548 : mapper_(&m), data_(&data)
552 Element
const &element()
const
562 FiniteElement
const &finiteElement()
const
564 return finiteElement_;
568 void bind(Element
const &e)
571 finiteElement_.bind((*data_)[mapper_->index(e)], *element_);
572 this->setSize(finiteElement_.size());
576 unsigned int order()
const
578 return finiteElement_.localBasis().order();
582 FiniteElement finiteElement_;
583 Element
const* element_;
584 Mapper
const* mapper_;
585 std::vector<std::bitset<3>>
const* data_;
598 template<
class GV,
class R>
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;
609 static constexpr auto morleyMapperLayout(Dune::GeometryType type,
int gridDim)
611 assert(gridDim == 2);
616 if ((type.isTriangle()) )
630 using Node = MorleyNode<GridView, R>;
634 :
Base(gv, morleyMapperLayout)
635 , mapper_({gv, mcmgElementLayout()})
637 data_ = Impl::computeEdgeOrientations(mapper_);
644 data_ = Impl::computeEdgeOrientations(mapper_);
652 return Node{mapper_, data_};
657 SubEntityMapper mapper_;
658 std::vector<std::bitset<3>> data_;
662 namespace BasisFactory
671 template<
class R =
double>
674 return [=](
auto const &gridView) {
675 return MorleyPreBasis<std::decay_t<
decltype(gridView)>, R>(gridView);
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