12#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_LOCAL_FINITE_ELEMENT_HH
13#define DUMUX_DISCRETIZATION_PQ1BUBBLE_LOCAL_FINITE_ELEMENT_HH
20#include <dune/common/fmatrix.hh>
21#include <dune/common/fvector.hh>
22#include <dune/common/exceptions.hh>
24#include <dune/geometry/type.hh>
25#include <dune/geometry/referenceelements.hh>
27#include <dune/localfunctions/common/localbasis.hh>
28#include <dune/localfunctions/common/localfiniteelementtraits.hh>
29#include <dune/localfunctions/common/localkey.hh>
31#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
32#include <dune/localfunctions/lagrange/lagrangecube.hh>
44template<
class D,
class R,
unsigned int dim, Dune::GeometryType::Id typeId>
47 using PQ1FiniteElement = std::conditional_t<
48 Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim),
49 Dune::LagrangeCubeLocalFiniteElement<D, R, dim, 1>,
50 Dune::LagrangeSimplexLocalFiniteElement<D, R, dim, 1>
52 static constexpr std::size_t numDofs
53 = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim) ? (1<<dim)+1 : (dim+1)+1;
55 using Traits = Dune::LocalBasisTraits<
56 D, dim, Dune::FieldVector<D, dim>,
57 R, 1, Dune::FieldVector<R, 1>,
58 Dune::FieldMatrix<R, 1, dim>
65 const auto& p1Basis = pq1FiniteElement_.localBasis();
66 const auto refElement = Dune::referenceElement<typename Traits::DomainFieldType, dim>(
type());
67 const auto&
center = refElement.position(0, 0);
68 p1Basis.evaluateFunction(
center, pq1AtCenter_);
74 static constexpr unsigned int size()
81 std::vector<typename Traits::RangeType>& out)
const
84 const auto& p1Basis = pq1FiniteElement_.localBasis();
85 p1Basis.evaluateFunction(x, out);
86 const auto bubble = evaluateBubble_(x);
89 for (
int i = 0; i < numDofs-1; ++i)
90 out[i] -= pq1AtCenter_[i]*out.back();
97 std::vector<typename Traits::JacobianType>& out)
const
100 const auto& p1Basis = pq1FiniteElement_.localBasis();
101 p1Basis.evaluateJacobian(x, out);
103 std::vector<typename Traits::RangeType> shapeValues;
104 p1Basis.evaluateFunction(x, shapeValues);
106 const auto bubbleJacobian = evaluateBubbleJacobian_(x);
108 for (
int i = 0; i < numDofs-1; ++i)
109 for (
int k = 0; k < dim; ++k)
110 out[i][0][k] -= pq1AtCenter_[i]*bubbleJacobian[0][k];
113 out.back() = bubbleJacobian;
123 const typename Traits::DomainType& in,
124 std::vector<typename Traits::RangeType>& out)
const
126 DUNE_THROW(Dune::NotImplemented,
"Partial derivatives");
133 static constexpr unsigned int order()
141 static constexpr Dune::GeometryType
type()
147 typename Traits::RangeType evaluateBubble_(
const typename Traits::DomainType& x)
const
149 if constexpr (
type() == Dune::GeometryTypes::simplex(dim))
151 if constexpr (dim == 2)
152 return 27*x[0]*x[1]*(1-x[0]-x[1]);
153 else if constexpr (dim == 3)
154 return 256*x[0]*x[1]*x[2]*(1-x[0]-x[1]-x[2]);
156 else if constexpr (
type() == Dune::GeometryTypes::cube(dim))
158 if constexpr (dim == 2)
159 return 16*x[0]*x[1]*(1-x[0])*(1-x[1]);
160 else if constexpr (dim == 3)
161 return 64*x[0]*x[1]*x[2]*(1-x[0])*(1-x[1])*(1-x[2]);
164 DUNE_THROW(Dune::NotImplemented,
"Bubble function for " <<
type());
168 typename Traits::JacobianType evaluateBubbleJacobian_(
const typename Traits::DomainType& x)
const
170 if constexpr (
type() == Dune::GeometryTypes::simplex(dim))
172 if constexpr (dim == 2)
173 return {{27*(x[1]*(1-x[0]-x[1]) - x[0]*x[1]),
174 27*(x[0]*(1-x[0]-x[1]) - x[0]*x[1])}};
175 else if constexpr (dim == 3)
176 return {{256*(x[1]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]),
177 256*(x[0]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]),
178 256*(x[0]*x[1]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2])}};
180 else if constexpr (
type() == Dune::GeometryTypes::cube(dim))
182 if constexpr (dim == 2)
183 return {{16*(x[1]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[1])),
184 16*(x[0]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[0]))}};
185 else if constexpr (dim == 3)
186 return {{64*(x[1]*x[2]*(1-x[0])*(1-x[1])*(1-x[2]) - x[0]*x[1]*x[2]*(1-x[1]))*(1-x[2]),
187 64*(x[0]*x[2]*(1-x[0])*(1-x[1])*(1-x[2]) - x[0]*x[1]*x[2]*(1-x[0]))*(1-x[2]),
188 64*(x[0]*x[1]*(1-x[0])*(1-x[1])*(1-x[2]) - x[0]*x[1]*x[2]*(1-x[0]))*(1-x[1])}};
191 DUNE_THROW(Dune::NotImplemented,
"Bubble function for " <<
type() <<
" dim = " << dim);
194 PQ1FiniteElement pq1FiniteElement_;
195 std::vector<typename Traits::RangeType> pq1AtCenter_;
203template<
int dim, Dune::GeometryType::Id typeId>
206 static constexpr std::size_t numDofs
207 = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim) ? (1<<dim)+1 : (dim+1)+1;
218 for (std::size_t i=0; i<
size()-1; i++)
219 localKeys_[i] = Dune::LocalKey(i, dim, 0);
222 localKeys_.back() = Dune::LocalKey(0, 0, 0);
226 static constexpr std::size_t
size ()
230 const Dune::LocalKey&
localKey (std::size_t i)
const
231 {
return localKeys_[i]; }
234 std::array<Dune::LocalKey, numDofs> localKeys_;
242template<
class LocalBasis>
254 template<
typename F,
typename C>
257 constexpr auto dim = LocalBasis::Traits::dimDomain;
259 out.resize(LocalBasis::size());
261 const auto refElement = Dune::referenceElement<typename LocalBasis::Traits::DomainFieldType,dim>(LocalBasis::type());
264 for (
int i = 0; i < refElement.size(dim); ++i)
265 out[i] = f(refElement.position(i, dim));
266 out.back() = f(refElement.position(0, 0));
282template<
class D,
class R,
int dim, Dune::GeometryType::Id typeId>
289 static constexpr Dune::GeometryType gt = Dune::GeometryType{ typeId };
291 gt == Dune::GeometryTypes::cube(dim) || gt == Dune::GeometryTypes::simplex(dim),
292 "Only implemented for cubes and simplices"
296 using Traits = Dune::LocalFiniteElementTraits<Basis, Coefficients, Interpolation>;
311 return coefficients_;
319 return interpolation_;
325 static constexpr std::size_t
size()
333 static constexpr Dune::GeometryType
type()
335 return Traits::LocalBasisType::type();
340 Coefficients coefficients_;
341 Interpolation interpolation_;
P1/Q1 + Bubble on the reference element.
Definition pq1bubblelocalfiniteelement.hh:46
void evaluateJacobian(const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const
Evaluate the Jacobians of all shape functions.
Definition pq1bubblelocalfiniteelement.hh:96
static constexpr Dune::GeometryType type()
The reference element type.
Definition pq1bubblelocalfiniteelement.hh:141
static constexpr unsigned int order()
Evaluate the Jacobians of all shape functions we are actually cubic/quartic but cannot represent all ...
Definition pq1bubblelocalfiniteelement.hh:133
void evaluateFunction(const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition pq1bubblelocalfiniteelement.hh:80
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of any order of all shape functions.
Definition pq1bubblelocalfiniteelement.hh:122
Dune::LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition pq1bubblelocalfiniteelement.hh:55
PQ1BubbleLocalBasis()
Definition pq1bubblelocalfiniteelement.hh:61
static constexpr unsigned int size()
Number of shape functions (one for each vertex and one in the element)
Definition pq1bubblelocalfiniteelement.hh:74
Associations of the P1/Q1 + Bubble degrees of freedom to the facets of the reference element.
Definition pq1bubblelocalfiniteelement.hh:205
static constexpr std::size_t size()
Number of coefficients.
Definition pq1bubblelocalfiniteelement.hh:226
PQ1BubbleLocalCoefficients()
Definition pq1bubblelocalfiniteelement.hh:210
const Dune::LocalKey & localKey(std::size_t i) const
Get i-th local key.
Definition pq1bubblelocalfiniteelement.hh:230
Evaluate the degrees of freedom of a P1 + Bubble basis.
Definition pq1bubblelocalfiniteelement.hh:244
void interpolate(const F &f, std::vector< C > &out) const
Evaluate a given function at the vertices and the cell midpoint.
Definition pq1bubblelocalfiniteelement.hh:255
P1/Q1 + Bubble finite element.
Definition pq1bubblelocalfiniteelement.hh:284
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition pq1bubblelocalfiniteelement.hh:317
Dune::LocalFiniteElementTraits< Basis, Coefficients, Interpolation > Traits
Definition pq1bubblelocalfiniteelement.hh:296
static constexpr std::size_t size()
The number of coefficients in the basis.
Definition pq1bubblelocalfiniteelement.hh:325
static constexpr Dune::GeometryType type()
The reference element type that the local finite element is defined on.
Definition pq1bubblelocalfiniteelement.hh:333
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition pq1bubblelocalfiniteelement.hh:301
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition pq1bubblelocalfiniteelement.hh:309
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition center.hh:24
Distance implementation details.
Definition cvfelocalresidual.hh:25