version 3.11-dev
pq1bubblelocalfiniteelement.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_LOCAL_FINITE_ELEMENT_HH
13#define DUMUX_DISCRETIZATION_PQ1BUBBLE_LOCAL_FINITE_ELEMENT_HH
14
15#include <array>
16#include <vector>
17#include <numeric>
18#include <algorithm>
19
20#include <dune/common/fmatrix.hh>
21#include <dune/common/fvector.hh>
22#include <dune/common/exceptions.hh>
23
24#include <dune/geometry/type.hh>
25#include <dune/geometry/referenceelements.hh>
26
27#include <dune/localfunctions/common/localbasis.hh>
28#include <dune/localfunctions/common/localfiniteelementtraits.hh>
29#include <dune/localfunctions/common/localkey.hh>
30
31#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
32#include <dune/localfunctions/lagrange/lagrangecube.hh>
33
34namespace Dumux::Detail {
35
45template<class D, class R, unsigned int dim, Dune::GeometryType::Id typeId, std::size_t numCubeBubbleDofs>
47{
48 using PQ1FiniteElement = std::conditional_t<
49 Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim),
50 Dune::LagrangeCubeLocalFiniteElement<D, R, dim, 1>,
51 Dune::LagrangeSimplexLocalFiniteElement<D, R, dim, 1>
52 >;
53 static constexpr bool isCube = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim);
54 static constexpr std::size_t numVertexDofs = isCube ? (1<<dim) : (dim+1);
55 static constexpr std::size_t numBubbleDofs = isCube ? numCubeBubbleDofs : 1;
56 static constexpr std::size_t numDofs = numVertexDofs + numBubbleDofs;
57public:
58 using Traits = Dune::LocalBasisTraits<
59 D, dim, Dune::FieldVector<D, dim>,
60 R, 1, Dune::FieldVector<R, 1>,
61 Dune::FieldMatrix<R, 1, dim>
62 >;
63
65 : pq1FiniteElement_{}
66 {
67 // precompute center values for normalization
68 const auto& p1Basis = pq1FiniteElement_.localBasis();
69 const auto refElement = Dune::referenceElement<typename Traits::DomainFieldType, dim>(type());
70 const auto& center = refElement.position(0, 0);
71 p1Basis.evaluateFunction(center, pq1AtCenter_);
72 }
73
77 static constexpr unsigned int size()
78 { return numDofs; }
79
83 void evaluateFunction(const typename Traits::DomainType& x,
84 std::vector<typename Traits::RangeType>& out) const
85 {
86 out.reserve(size());
87 const auto& p1Basis = pq1FiniteElement_.localBasis();
88 p1Basis.evaluateFunction(x, out);
89 const auto bubble = evaluateBubble_(x);
90
91 out.resize(size());
92 out[numVertexDofs] = bubble;
93
94 if constexpr (numBubbleDofs == 2)
95 out[numVertexDofs+1] = scaling_(x)*bubble;
96
97 for (std::size_t i = 0; i < numVertexDofs; ++i)
98 out[i] -= pq1AtCenter_[i]*out[numVertexDofs];
99 }
100
104 void evaluateJacobian(const typename Traits::DomainType& x,
105 std::vector<typename Traits::JacobianType>& out) const
106 {
107 out.reserve(size());
108 const auto& p1Basis = pq1FiniteElement_.localBasis();
109 p1Basis.evaluateJacobian(x, out);
110
111 std::vector<typename Traits::RangeType> shapeValues;
112 p1Basis.evaluateFunction(x, shapeValues);
113
114 const auto bubble = evaluateBubble_(x);
115 const auto bubbleJacobian = evaluateBubbleJacobian_(x);
116
117 for (std::size_t i = 0; i < numVertexDofs; ++i)
118 for (int k = 0; k < dim; ++k)
119 out[i][0][k] -= pq1AtCenter_[i]*bubbleJacobian[0][k];
120
121 out.resize(size());
122 out[numVertexDofs] = bubbleJacobian;
123 if constexpr (numBubbleDofs == 2)
124 out[numVertexDofs+1] = scaling_(x)*bubbleJacobian + bubble[0]*scalingJacobian_(x);
125 }
126
133 void partial(const std::array<unsigned int, dim>& order,
134 const typename Traits::DomainType& in,
135 std::vector<typename Traits::RangeType>& out) const
136 {
137 DUNE_THROW(Dune::NotImplemented, "Partial derivatives");
138 }
139
144 static constexpr unsigned int order()
145 {
146 return 1;
147 }
148
152 static constexpr Dune::GeometryType type()
153 {
154 return { typeId };
155 }
156private:
157 // evaluate bubble function at x
158 typename Traits::RangeType evaluateBubble_(const typename Traits::DomainType& x) const
159 {
160 if constexpr (type() == Dune::GeometryTypes::simplex(dim))
161 {
162 if constexpr (dim == 2)
163 return 27*x[0]*x[1]*(1-x[0]-x[1]);
164 else if constexpr (dim == 3)
165 return 256*x[0]*x[1]*x[2]*(1-x[0]-x[1]-x[2]);
166 }
167 else if constexpr (type() == Dune::GeometryTypes::cube(dim))
168 {
169 if constexpr (dim == 2)
170 return 16*x[0]*x[1]*(1-x[0])*(1-x[1]);
171 else if constexpr (dim == 3)
172 return 64*x[0]*x[1]*x[2]*(1-x[0])*(1-x[1])*(1-x[2]);
173 }
174 else
175 DUNE_THROW(Dune::NotImplemented, "Bubble function for " << type());
176 }
177
178 // evaluate bubble function at x
179 typename Traits::JacobianType evaluateBubbleJacobian_(const typename Traits::DomainType& x) const
180 {
181 if constexpr (type() == Dune::GeometryTypes::simplex(dim))
182 {
183 if constexpr (dim == 2)
184 return {{27*(x[1]*(1-x[0]-x[1]) - x[0]*x[1]),
185 27*(x[0]*(1-x[0]-x[1]) - x[0]*x[1])}};
186 else if constexpr (dim == 3)
187 return {{256*(x[1]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]),
188 256*(x[0]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]),
189 256*(x[0]*x[1]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2])}};
190 }
191 else if constexpr (type() == Dune::GeometryTypes::cube(dim))
192 {
193 if constexpr (dim == 2)
194 return {{16*(x[1]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[1])),
195 16*(x[0]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[0]))}};
196 else if constexpr (dim == 3)
197 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]),
198 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]),
199 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])}};
200 }
201 else
202 DUNE_THROW(Dune::NotImplemented, "Bubble function for " << type() << " dim = " << dim);
203 }
204
205 auto scaling_(const typename Traits::DomainType& x) const
206 {
207 if constexpr (dim == 2)
208 return (x[0]-0.5) + (x[1]-0.5);
209 else if constexpr (dim == 3)
210 return (x[0]-0.5) + (x[1]-0.5) + (x[2]-0.5);
211 else
212 DUNE_THROW(Dune::NotImplemented, "Bubble function function for dim = " << dim);
213 }
214
215 typename Traits::JacobianType scalingJacobian_(const typename Traits::DomainType& x) const
216 {
217 if constexpr (dim == 2)
218 return {{1.0, 1.0}};
219 else if constexpr (dim == 3)
220 return {{1.0, 1.0, 1.0}};
221 else
222 DUNE_THROW(Dune::NotImplemented, "Bubble function function for dim = " << dim);
223 }
224
225 PQ1FiniteElement pq1FiniteElement_;
226 std::vector<typename Traits::RangeType> pq1AtCenter_;
227};
228
235template<int dim, Dune::GeometryType::Id typeId, std::size_t numCubeBubbleDofs>
237{
238 static constexpr bool isCube = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim);
239 static constexpr std::size_t numVertexDofs = isCube ? (1<<dim) : (dim+1);
240 static constexpr std::size_t numBubbleDofs = isCube ? numCubeBubbleDofs : 1;
241 static constexpr std::size_t numDofs = numVertexDofs + numBubbleDofs;
242public:
243
245 {
246 // Dune::LocalKey is constructed with
247 // - the local index with respect to the subentity type
248 // - the codim of the subentity the dof is attached to
249 // - the local number of the function to evaluate
250
251 // vertex dofs
252 for (std::size_t i = 0; i < numVertexDofs; ++i)
253 localKeys_[i] = Dune::LocalKey(i, dim, 0);
254
255 // cell dof(s)
256 for (std::size_t i = 0; i < numBubbleDofs; ++i)
257 localKeys_[numVertexDofs + i] = Dune::LocalKey(0, 0, i);
258 }
259
261 static constexpr std::size_t size ()
262 { return numDofs; }
263
265 const Dune::LocalKey& localKey (std::size_t i) const
266 { return localKeys_[i]; }
267
268private:
269 std::array<Dune::LocalKey, numDofs> localKeys_;
270};
271
277template<class LocalBasis>
279{
280public:
289 template<typename F, typename C>
290 void interpolate (const F& f, std::vector<C>& out) const
291 {
292 constexpr auto dim = LocalBasis::Traits::dimDomain;
293 if constexpr (LocalBasis::type() == Dune::GeometryTypes::cube(dim))
294 DUNE_THROW(Dune::NotImplemented, "Interpolation for Q1 bubble not implemented yet");
295
296 out.resize(LocalBasis::size());
297
298 const auto refElement = Dune::referenceElement<typename LocalBasis::Traits::DomainFieldType,dim>(LocalBasis::type());
299
300 // Evaluate at the vertices and at the center
301 for (int i = 0; i < refElement.size(dim); ++i)
302 out[i] = f(refElement.position(i, dim));
303 out.back() = f(refElement.position(0, 0));
304 }
305};
306
307} // end namespace Dumux::Detail
308
309namespace Dumux {
310
320template<class D, class R, int dim, Dune::GeometryType::Id typeId, std::size_t numCubeBubbleDofs = 1>
322{
323 static_assert(numCubeBubbleDofs == 1 || numCubeBubbleDofs == 2,
324 "P1/Q1 bubble FE supports numCubeBubbleDofs = 1 or 2");
328
329 static constexpr Dune::GeometryType gt = Dune::GeometryType{ typeId };
330 static_assert(
331 gt == Dune::GeometryTypes::cube(dim) || gt == Dune::GeometryTypes::simplex(dim),
332 "Only implemented for cubes and simplices"
333 );
334
335public:
336 using Traits = Dune::LocalFiniteElementTraits<Basis, Coefficients, Interpolation>;
337
341 const typename Traits::LocalBasisType& localBasis() const
342 {
343 return basis_;
344 }
345
349 const typename Traits::LocalCoefficientsType& localCoefficients() const
350 {
351 return coefficients_;
352 }
353
357 const typename Traits::LocalInterpolationType& localInterpolation() const
358 {
359 return interpolation_;
360 }
361
365 static constexpr std::size_t size()
366 {
367 return Basis::size();
368 }
369
373 static constexpr Dune::GeometryType type()
374 {
375 return Traits::LocalBasisType::type();
376 }
377
378private:
379 Basis basis_;
380 Coefficients coefficients_;
381 Interpolation interpolation_;
382};
383
384} // end namespace Dumux
385
386#endif
P1/Q1 + Bubble functions on the reference element.
Definition: pq1bubblelocalfiniteelement.hh:47
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:133
static constexpr unsigned int order()
Evaluate the Jacobians of all shape functions we are actually cubic/quartic but cannot represent all ...
Definition: pq1bubblelocalfiniteelement.hh:144
void evaluateJacobian(const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const
Evaluate the Jacobians of all shape functions.
Definition: pq1bubblelocalfiniteelement.hh:104
static constexpr Dune::GeometryType type()
The reference element type.
Definition: pq1bubblelocalfiniteelement.hh:152
static constexpr unsigned int size()
Number of shape functions (one for each vertex and one in the element)
Definition: pq1bubblelocalfiniteelement.hh:77
void evaluateFunction(const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: pq1bubblelocalfiniteelement.hh:83
PQ1BubbleLocalBasis()
Definition: pq1bubblelocalfiniteelement.hh:64
Dune::LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition: pq1bubblelocalfiniteelement.hh:62
Associations of the P1/Q1 + Bubble degrees of freedom to the facets of the reference element.
Definition: pq1bubblelocalfiniteelement.hh:237
const Dune::LocalKey & localKey(std::size_t i) const
Get i-th local key.
Definition: pq1bubblelocalfiniteelement.hh:265
PQ1BubbleLocalCoefficients()
Definition: pq1bubblelocalfiniteelement.hh:244
static constexpr std::size_t size()
Number of coefficients.
Definition: pq1bubblelocalfiniteelement.hh:261
Evaluate the degrees of freedom of a P1 + Bubble basis.
Definition: pq1bubblelocalfiniteelement.hh:279
void interpolate(const F &f, std::vector< C > &out) const
Evaluate a given function at the vertices and the cell midpoint.
Definition: pq1bubblelocalfiniteelement.hh:290
P1/Q1 + Bubble finite element.
Definition: pq1bubblelocalfiniteelement.hh:322
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: pq1bubblelocalfiniteelement.hh:357
static constexpr std::size_t size()
The number of coefficients in the basis.
Definition: pq1bubblelocalfiniteelement.hh:365
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: pq1bubblelocalfiniteelement.hh:341
static constexpr Dune::GeometryType type()
The reference element type that the local finite element is defined on.
Definition: pq1bubblelocalfiniteelement.hh:373
Dune::LocalFiniteElementTraits< Basis, Coefficients, Interpolation > Traits
Definition: pq1bubblelocalfiniteelement.hh:336
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: pq1bubblelocalfiniteelement.hh:349
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Definition: cvfelocalresidual.hh:30
Definition: adapt.hh:17