Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup PQ1BubbleDiscretization | ||
10 | * \brief Evaluate P1/Q1 basis with bubble function | ||
11 | */ | ||
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 | |||
35 | namespace Dumux::Detail { | ||
36 | |||
37 | /*! | ||
38 | * \brief P1/Q1 + Bubble on the reference element | ||
39 | * | ||
40 | * \tparam D Type to represent the field in the domain | ||
41 | * \tparam R Type to represent the field in the range | ||
42 | * \tparam dim Dimension of the domain element | ||
43 | * \tparam typeId The geometry type | ||
44 | */ | ||
45 | template<class D, class R, unsigned int dim, Dune::GeometryType::Id typeId> | ||
46 | class PQ1BubbleLocalBasis | ||
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; | ||
55 | public: | ||
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 | |||
62 | 52 | PQ1BubbleLocalBasis() | |
63 |
1/4✓ Branch 2 taken 26 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
52 | : pq1FiniteElement_{} |
64 | { | ||
65 | // precompute center values for normalization | ||
66 |
1/2✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
|
52 | const auto& p1Basis = pq1FiniteElement_.localBasis(); |
67 |
1/4✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
52 | const auto refElement = Dune::referenceElement<typename Traits::DomainFieldType, dim>(type()); |
68 |
1/2✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
|
52 | const auto& center = refElement.position(0, 0); |
69 |
1/2✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
|
52 | p1Basis.evaluateFunction(center, pq1AtCenter_); |
70 | 52 | } | |
71 | |||
72 | /*! | ||
73 | * \brief Number of shape functions (one for each vertex and one in the element) | ||
74 | */ | ||
75 | static constexpr unsigned int size() | ||
76 | { return numDofs; } | ||
77 | |||
78 | /*! | ||
79 | * \brief Evaluate all shape functions | ||
80 | */ | ||
81 | 33966068 | void evaluateFunction(const typename Traits::DomainType& x, | |
82 | std::vector<typename Traits::RangeType>& out) const | ||
83 | { | ||
84 | 33966068 | out.reserve(size()); | |
85 | 33966068 | const auto& p1Basis = pq1FiniteElement_.localBasis(); | |
86 | 33966068 | p1Basis.evaluateFunction(x, out); | |
87 | 67932136 | const auto bubble = evaluateBubble_(x); | |
88 | 33966068 | out.resize(size()); | |
89 | 33966068 | out.back() = bubble; | |
90 |
2/2✓ Branch 0 taken 59406428 times.
✓ Branch 1 taken 16983034 times.
|
152778924 | for (int i = 0; i < numDofs-1; ++i) |
91 | 712877136 | out[i] -= pq1AtCenter_[i]*out.back(); | |
92 | 33966068 | } | |
93 | |||
94 | /*! | ||
95 | * \brief Evaluate the Jacobians of all shape functions | ||
96 | */ | ||
97 | 1755172 | void evaluateJacobian(const typename Traits::DomainType& x, | |
98 | std::vector<typename Traits::JacobianType>& out) const | ||
99 | { | ||
100 | 1755172 | out.reserve(size()); | |
101 | 1755172 | const auto& p1Basis = pq1FiniteElement_.localBasis(); | |
102 | 1755172 | p1Basis.evaluateJacobian(x, out); | |
103 | |||
104 |
1/4✓ Branch 1 taken 877586 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
1755172 | std::vector<typename Traits::RangeType> shapeValues; |
105 |
1/2✓ Branch 1 taken 877586 times.
✗ Branch 2 not taken.
|
1755172 | p1Basis.evaluateFunction(x, shapeValues); |
106 | |||
107 | 1755172 | const auto bubbleJacobian = evaluateBubbleJacobian_(x); | |
108 | |||
109 |
2/2✓ Branch 1 taken 3267252 times.
✓ Branch 2 taken 877586 times.
|
8289676 | for (int i = 0; i < numDofs-1; ++i) |
110 |
2/2✓ Branch 0 taken 6741152 times.
✓ Branch 1 taken 3267252 times.
|
20016808 | for (int k = 0; k < dim; ++k) |
111 | 107858432 | out[i][0][k] -= pq1AtCenter_[i]*bubbleJacobian[0][k]; | |
112 | |||
113 |
1/2✓ Branch 1 taken 877586 times.
✗ Branch 2 not taken.
|
1755172 | out.resize(size()); |
114 |
2/4✓ Branch 0 taken 877586 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 877586 times.
✗ Branch 3 not taken.
|
3510344 | out.back() = bubbleJacobian; |
115 | 1755172 | } | |
116 | |||
117 | /** \brief Evaluate partial derivatives of any order of all shape functions | ||
118 | * | ||
119 | * \param order Order of the partial derivatives, in the classic multi-index notation | ||
120 | * \param in Position where to evaluate the derivatives | ||
121 | * \param[out] out The desired partial derivatives | ||
122 | */ | ||
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 | |||
130 | /*! | ||
131 | * \brief Evaluate the Jacobians of all shape functions | ||
132 | * we are actually cubic/quartic but cannot represent all cubic/quartic polynomials | ||
133 | */ | ||
134 | static constexpr unsigned int order() | ||
135 | { | ||
136 | return 1; | ||
137 | } | ||
138 | |||
139 | /*! | ||
140 | * \brief The reference element type | ||
141 | */ | ||
142 | static constexpr Dune::GeometryType type() | ||
143 | { | ||
144 | return { typeId }; | ||
145 | } | ||
146 | private: | ||
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 | 42628540 | return 27*x[0]*x[1]*(1-x[0]-x[1]); | |
154 | else if constexpr (dim == 3) | ||
155 | 5278434 | 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 | 38516320 | 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 | 1755172 | 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 | 2430920 | return {{27*(x[1]*(1-x[0]-x[1]) - x[0]*x[1]), | |
175 | 3403288 | 27*(x[0]*(1-x[0]-x[1]) - x[0]*x[1])}}; | |
176 | else if constexpr (dim == 3) | ||
177 | 826592 | return {{256*(x[1]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]), | |
178 | 723268 | 256*(x[0]*x[2]*(1-x[0]-x[1]-x[2]) - x[0]*x[1]*x[2]), | |
179 | 1033240 | 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 | 6993984 | return {{16*(x[1]*(1-x[0])*(1-x[1]) - x[0]*x[1]*(1-x[1])), | |
185 | 9325312 | 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 | |||
199 | /*! | ||
200 | * \brief Associations of the P1/Q1 + Bubble degrees of freedom to the facets of the reference element | ||
201 | * \tparam dim Dimension of the reference element | ||
202 | * \tparam typeId The geometry type | ||
203 | */ | ||
204 | template<int dim, Dune::GeometryType::Id typeId> | ||
205 | class PQ1BubbleLocalCoefficients | ||
206 | { | ||
207 | static constexpr std::size_t numDofs | ||
208 | = Dune::GeometryType{ typeId } == Dune::GeometryTypes::cube(dim) ? (1<<dim)+1 : (dim+1)+1; | ||
209 | public: | ||
210 | |||
211 | PQ1BubbleLocalCoefficients() | ||
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 | |||
226 | //! Number of coefficients | ||
227 | static constexpr std::size_t size () | ||
228 | { return numDofs; } | ||
229 | |||
230 | //! Get i-th local key | ||
231 | const Dune::LocalKey& localKey (std::size_t i) const | ||
232 | 8501824 | { return localKeys_[i]; } | |
233 | |||
234 | private: | ||
235 | std::array<Dune::LocalKey, numDofs> localKeys_; | ||
236 | }; | ||
237 | |||
238 | /*! | ||
239 | * \brief Evaluate the degrees of freedom of a P1 + Bubble basis | ||
240 | * | ||
241 | * \tparam LocalBasis The corresponding set of shape functions | ||
242 | */ | ||
243 | template<class LocalBasis> | ||
244 | class PQ1BubbleLocalInterpolation | ||
245 | { | ||
246 | public: | ||
247 | /*! | ||
248 | * \brief Evaluate a given function at the vertices and the cell midpoint | ||
249 | * | ||
250 | * \tparam F Type of function to evaluate | ||
251 | * \tparam C Type used for the values of the function | ||
252 | * \param[in] f Function to evaluate (call operator that gets a local position) | ||
253 | * \param[out] out Array of function values | ||
254 | */ | ||
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 | |||
273 | namespace Dumux { | ||
274 | |||
275 | /*! | ||
276 | * \brief P1/Q1 + Bubble finite element | ||
277 | * | ||
278 | * \tparam D type used for domain coordinates | ||
279 | * \tparam R type used for function values | ||
280 | * \tparam dim dimension of the reference element | ||
281 | * \tparam typeId The geometry type | ||
282 | */ | ||
283 | template<class D, class R, int dim, Dune::GeometryType::Id typeId> | ||
284 |
2/4✓ Branch 1 taken 13 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 13 times.
✗ Branch 5 not taken.
|
78 | class PQ1BubbleLocalFiniteElement |
285 | { | ||
286 | using Basis = Detail::PQ1BubbleLocalBasis<D, R, dim, typeId>; | ||
287 | using Coefficients = Detail::PQ1BubbleLocalCoefficients<dim, typeId>; | ||
288 | using Interpolation = Detail::PQ1BubbleLocalInterpolation<Basis>; | ||
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 | |||
296 | public: | ||
297 | using Traits = Dune::LocalFiniteElementTraits<Basis, Coefficients, Interpolation>; | ||
298 | |||
299 | /*! | ||
300 | * \brief Returns the local basis, i.e., the set of shape functions | ||
301 | */ | ||
302 | const typename Traits::LocalBasisType& localBasis() const | ||
303 | { | ||
304 | 26 | return basis_; | |
305 | } | ||
306 | |||
307 | /*! | ||
308 | * \brief Returns the assignment of the degrees of freedom to the element subentities | ||
309 | */ | ||
310 | const typename Traits::LocalCoefficientsType& localCoefficients() const | ||
311 | { | ||
312 | 26 | return coefficients_; | |
313 | } | ||
314 | |||
315 | /*! | ||
316 | * \brief Returns object that evaluates degrees of freedom | ||
317 | */ | ||
318 | const typename Traits::LocalInterpolationType& localInterpolation() const | ||
319 | { | ||
320 | 26 | return interpolation_; | |
321 | } | ||
322 | |||
323 | /*! | ||
324 | * \brief The number of coefficients in the basis | ||
325 | */ | ||
326 | static constexpr std::size_t size() | ||
327 | { | ||
328 | return Basis::size(); | ||
329 | } | ||
330 | |||
331 | /*! | ||
332 | * \brief The reference element type that the local finite element is defined on | ||
333 | */ | ||
334 | static constexpr Dune::GeometryType type() | ||
335 | { | ||
336 | return Traits::LocalBasisType::type(); | ||
337 | } | ||
338 | |||
339 | private: | ||
340 | Basis basis_; | ||
341 | Coefficients coefficients_; | ||
342 | Interpolation interpolation_; | ||
343 | }; | ||
344 | |||
345 | } // end namespace Dumux | ||
346 | |||
347 | #endif | ||
348 |