version 3.8
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-FileCopyrightInfo: 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/localinterpolation.hh>
30#include <dune/localfunctions/common/localkey.hh>
31
32#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
33#include <dune/localfunctions/lagrange/lagrangecube.hh>
34
35namespace Dumux::Detail {
36
45template<class D, class R, unsigned int dim, Dune::GeometryType::Id typeId>
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 std::size_t numDofs
54 = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim) ? (1<<dim)+1 : (dim+1)+1;
55public:
56 using Traits = Dune::LocalBasisTraits<
57 D, dim, Dune::FieldVector<D, dim>,
58 R, 1, Dune::FieldVector<R, 1>,
59 Dune::FieldMatrix<R, 1, dim>
60 >;
61
63 : pq1FiniteElement_{}
64 {
65 // precompute center values for normalization
66 const auto& p1Basis = pq1FiniteElement_.localBasis();
67 const auto refElement = Dune::referenceElement<typename Traits::DomainFieldType, dim>(type());
68 const auto& center = refElement.position(0, 0);
69 p1Basis.evaluateFunction(center, pq1AtCenter_);
70 }
71
75 static constexpr unsigned int size()
76 { return numDofs; }
77
81 void evaluateFunction(const typename Traits::DomainType& x,
82 std::vector<typename Traits::RangeType>& out) const
83 {
84 out.reserve(size());
85 const auto& p1Basis = pq1FiniteElement_.localBasis();
86 p1Basis.evaluateFunction(x, out);
87 const auto bubble = evaluateBubble_(x);
88 out.resize(size());
89 out.back() = bubble;
90 for (int i = 0; i < numDofs-1; ++i)
91 out[i] -= pq1AtCenter_[i]*out.back();
92 }
93
97 void evaluateJacobian(const typename Traits::DomainType& x,
98 std::vector<typename Traits::JacobianType>& out) const
99 {
100 out.reserve(size());
101 const auto& p1Basis = pq1FiniteElement_.localBasis();
102 p1Basis.evaluateJacobian(x, out);
103
104 std::vector<typename Traits::RangeType> shapeValues;
105 p1Basis.evaluateFunction(x, shapeValues);
106
107 const auto bubbleJacobian = evaluateBubbleJacobian_(x);
108
109 for (int i = 0; i < numDofs-1; ++i)
110 for (int k = 0; k < dim; ++k)
111 out[i][0][k] -= pq1AtCenter_[i]*bubbleJacobian[0][k];
112
113 out.resize(size());
114 out.back() = bubbleJacobian;
115 }
116
123 void partial(const std::array<unsigned int, dim>& order,
124 const typename Traits::DomainType& in,
125 std::vector<typename Traits::RangeType>& out) const
126 {
127 DUNE_THROW(Dune::NotImplemented, "Partial derivatives");
128 }
129
134 static constexpr unsigned int order()
135 {
136 return 1;
137 }
138
142 static constexpr Dune::GeometryType type()
143 {
144 return { typeId };
145 }
146private:
147 // evaluate bubble function at x
148 typename Traits::RangeType evaluateBubble_(const typename Traits::DomainType& x) const
149 {
150 if constexpr (type() == Dune::GeometryTypes::simplex(dim))
151 {
152 if constexpr (dim == 2)
153 return 27*x[0]*x[1]*(1-x[0]-x[1]);
154 else if constexpr (dim == 3)
155 return 256*x[0]*x[1]*x[2]*(1-x[0]-x[1]-x[2]);
156 }
157 else if constexpr (type() == Dune::GeometryTypes::cube(dim))
158 {
159 if constexpr (dim == 2)
160 return 16*x[0]*x[1]*(1-x[0])*(1-x[1]);
161 else if constexpr (dim == 3)
162 return 64*x[0]*x[1]*x[2]*(1-x[0])*(1-x[1])*(1-x[2]);
163 }
164 else
165 DUNE_THROW(Dune::NotImplemented, "Bubble function for " << type());
166 }
167
168 // evaluate bubble function at x
169 typename Traits::JacobianType evaluateBubbleJacobian_(const typename Traits::DomainType& x) const
170 {
171 if constexpr (type() == Dune::GeometryTypes::simplex(dim))
172 {
173 if constexpr (dim == 2)
174 return {{27*(x[1]*(1-x[0]-x[1]) - x[0]*x[1]),
175 27*(x[0]*(1-x[0]-x[1]) - x[0]*x[1])}};
176 else if constexpr (dim == 3)
177 return {{256*(x[1]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]),
178 256*(x[0]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]),
179 256*(x[0]*x[1]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2])}};
180 }
181 else if constexpr (type() == Dune::GeometryTypes::cube(dim))
182 {
183 if constexpr (dim == 2)
184 return {{16*(x[1]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[1])),
185 16*(x[0]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[0]))}};
186 else if constexpr (dim == 3)
187 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]),
188 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]),
189 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])}};
190 }
191 else
192 DUNE_THROW(Dune::NotImplemented, "Bubble function for " << type() << " dim = " << dim);
193 }
194
195 PQ1FiniteElement pq1FiniteElement_;
196 std::vector<typename Traits::RangeType> pq1AtCenter_;
197};
198
204template<int dim, Dune::GeometryType::Id typeId>
206{
207 static constexpr std::size_t numDofs
208 = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim) ? (1<<dim)+1 : (dim+1)+1;
209public:
210
212 {
213 // Dune::LocalKey is constructed with
214 // - the local index with respect to the subentity type
215 // - the codim of the subentity the dof is attached to
216 // - the local number of the function to evaluate
217
218 // vertex dofs
219 for (std::size_t i=0; i<size()-1; i++)
220 localKeys_[i] = Dune::LocalKey(i, dim, 0);
221
222 // cell dof
223 localKeys_.back() = Dune::LocalKey(0, 0, 0);
224 }
225
227 static constexpr std::size_t size ()
228 { return numDofs; }
229
231 const Dune::LocalKey& localKey (std::size_t i) const
232 { return localKeys_[i]; }
233
234private:
235 std::array<Dune::LocalKey, numDofs> localKeys_;
236};
237
243template<class LocalBasis>
245{
246public:
255 template<typename F, typename C>
256 void interpolate (const F& f, std::vector<C>& out) const
257 {
258 constexpr auto dim = LocalBasis::Traits::dimDomain;
259
260 out.resize(LocalBasis::size());
261
262 const auto refElement = Dune::referenceElement<typename LocalBasis::Traits::DomainFieldType,dim>(LocalBasis::type());
263
264 // Evaluate at the vertices and at the center
265 for (int i = 0; i < refElement.size(dim); ++i)
266 out[i] = f(refElement.position(i, dim));
267 out.back() = f(refElement.position(0, 0));
268 }
269};
270
271} // end namespace Dumux::Detail
272
273namespace Dumux {
274
283template<class D, class R, int dim, Dune::GeometryType::Id typeId>
285{
289
290 static constexpr Dune::GeometryType gt = Dune::GeometryType{ typeId };
291 static_assert(
292 gt == Dune::GeometryTypes::cube(dim) || gt == Dune::GeometryTypes::simplex(dim),
293 "Only implemented for cubes and simplices"
294 );
295
296public:
297 using Traits = Dune::LocalFiniteElementTraits<Basis, Coefficients, Interpolation>;
298
302 const typename Traits::LocalBasisType& localBasis() const
303 {
304 return basis_;
305 }
306
310 const typename Traits::LocalCoefficientsType& localCoefficients() const
311 {
312 return coefficients_;
313 }
314
318 const typename Traits::LocalInterpolationType& localInterpolation() const
319 {
320 return interpolation_;
321 }
322
326 static constexpr std::size_t size()
327 {
328 return Basis::size();
329 }
330
334 static constexpr Dune::GeometryType type()
335 {
336 return Traits::LocalBasisType::type();
337 }
338
339private:
340 Basis basis_;
341 Coefficients coefficients_;
342 Interpolation interpolation_;
343};
344
345} // end namespace Dumux
346
347#endif
P1/Q1 + Bubble on the reference element.
Definition: pq1bubblelocalfiniteelement.hh:47
void evaluateJacobian(const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const
Evaluate the Jacobians of all shape functions.
Definition: pq1bubblelocalfiniteelement.hh:97
static constexpr Dune::GeometryType type()
The reference element type.
Definition: pq1bubblelocalfiniteelement.hh:142
Dune::LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition: pq1bubblelocalfiniteelement.hh:60
static constexpr unsigned int order()
Evaluate the Jacobians of all shape functions we are actually cubic/quartic but cannot represent all ...
Definition: pq1bubblelocalfiniteelement.hh:134
void evaluateFunction(const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: pq1bubblelocalfiniteelement.hh:81
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:123
PQ1BubbleLocalBasis()
Definition: pq1bubblelocalfiniteelement.hh:62
static constexpr unsigned int size()
Number of shape functions (one for each vertex and one in the element)
Definition: pq1bubblelocalfiniteelement.hh:75
Associations of the P1/Q1 + Bubble degrees of freedom to the facets of the reference element.
Definition: pq1bubblelocalfiniteelement.hh:206
static constexpr std::size_t size()
Number of coefficients.
Definition: pq1bubblelocalfiniteelement.hh:227
PQ1BubbleLocalCoefficients()
Definition: pq1bubblelocalfiniteelement.hh:211
const Dune::LocalKey & localKey(std::size_t i) const
Get i-th local key.
Definition: pq1bubblelocalfiniteelement.hh:231
Evaluate the degrees of freedom of a P1 + Bubble basis.
Definition: pq1bubblelocalfiniteelement.hh:245
void interpolate(const F &f, std::vector< C > &out) const
Evaluate a given function at the vertices and the cell midpoint.
Definition: pq1bubblelocalfiniteelement.hh:256
P1/Q1 + Bubble finite element.
Definition: pq1bubblelocalfiniteelement.hh:285
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: pq1bubblelocalfiniteelement.hh:318
Dune::LocalFiniteElementTraits< Basis, Coefficients, Interpolation > Traits
Definition: pq1bubblelocalfiniteelement.hh:297
static constexpr std::size_t size()
The number of coefficients in the basis.
Definition: pq1bubblelocalfiniteelement.hh:326
static constexpr Dune::GeometryType type()
The reference element type that the local finite element is defined on.
Definition: pq1bubblelocalfiniteelement.hh:334
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: pq1bubblelocalfiniteelement.hh:302
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: pq1bubblelocalfiniteelement.hh:310
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
Definition: adapt.hh:17