MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_pyramid.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_FE_PYRAMID
13#define MFEM_FE_PYRAMID
14
15#include "fe_base.hpp"
16
17namespace mfem
18{
19
20/** Base class for arbitrary order basis functions on pyramid-shaped elements
21
22 This base class provides a common class to store temporary vectors,
23 matrices, and tensors computed by various functions defined on
24 pyramid-shaped elements.
25
26 The function names defined here are chosen to reflect, as closely as
27 possible, those used in the paper "Orientation embedded high order
28 shape functions for the exact sequence elements of all shapes" by
29 Federico Fuentes, Brendan Keith, Leszek Demkowicz, and Sriram
30 Nagaraj, see https://doi.org/10.1016/j.camwa.2015.04.027.
31
32 @note Many of the functions below, e.g. lam1, lam2, etc. and related
33 functions, are singular or multi-valued at the apex of the pyramid. The
34 values returned near the apex are computed in the limit z->1 using
35 (x, y, z) = ((1-z)/2, (1-z)/2, z) i.e. along the line from the center
36 of the base to the apex.
37*/
39{
40private:
41#ifndef MFEM_THREAD_SAFE
42 mutable DenseMatrix phi_E_mtmp;
43 mutable Vector phi_Q_vtmp1;
44 mutable Vector phi_Q_vtmp2;
45 mutable DenseMatrix phi_Q_mtmp1;
46 mutable DenseMatrix phi_Q_mtmp2;
47 mutable Vector phi_T_vtmp1;
48 mutable Vector phi_T_vtmp2;
49 mutable Vector phi_T_vtmp3;
50 mutable Vector phi_T_vtmp4;
51 mutable DenseMatrix phi_T_mtmp1;
52 mutable Vector E_E_vtmp;
53 mutable Vector E_Q_vtmp;
54 mutable DenseMatrix E_Q_mtmp1;
55 mutable DenseMatrix E_Q_mtmp2;
56 mutable DenseMatrix E_Q_mtmp3;
57 mutable Vector E_T_vtmp1;
58 mutable Vector E_T_vtmp2;
59 mutable Vector E_T_vtmp3;
60 mutable DenseMatrix E_T_mtmp1;
61 mutable DenseMatrix E_T_mtmp2;
62 mutable DenseMatrix V_Q_mtmp1;
63 mutable DenseMatrix V_Q_mtmp2;
64 mutable Vector V_T_vtmp1;
65 mutable Vector V_T_vtmp2;
66 mutable DenseMatrix V_T_mtmp1;
67 mutable Vector VT_T_vtmp1;
68 mutable Vector VT_T_vtmp2;
69 mutable DenseMatrix VT_T_mtmp1;
70 mutable DenseTensor VT_T_ttmp1;
71 mutable Vector V_L_vtmp1;
72 mutable Vector V_L_vtmp2;
73 mutable DenseMatrix V_L_mtmp1;
74 mutable DenseMatrix V_L_mtmp2;
75 mutable Vector V_R_vtmp;
76 mutable DenseMatrix V_R_mtmp;
77#endif
78
79protected:
80 static constexpr real_t one = 1.0;
81 static constexpr real_t zero = 0.0;
82 static constexpr real_t apex_tol = 1e-8;
83
84public:
85 FuentesPyramid() = default;
86
87 static bool CheckZ(real_t z) { return std::abs(z - 1.0) > apex_tol; }
88
89 /// Pyramid "Affine" Coordinates
90 static real_t lam1(real_t x, real_t y, real_t z)
91 { return CheckZ(z) ? (1.0 - x - z) * (1.0 - y - z) / (1.0 - z): 0.0; }
92 static real_t lam2(real_t x, real_t y, real_t z)
93 { return CheckZ(z) ? x * (1.0 - y - z) / (1.0 - z): 0.0; }
94 static real_t lam3(real_t x, real_t y, real_t z)
95 { return CheckZ(z) ? x * y / (1.0 - z): 0.0; }
96 static real_t lam4(real_t x, real_t y, real_t z)
97 { return CheckZ(z) ? (1.0 - x - z) * y / (1.0 - z): 0.0; }
98 static real_t lam5(real_t x, real_t y, real_t z)
99 { return CheckZ(z) ? z : 1.0; }
100
101 /// Gradients of the "Affine" Coordinates
102 static Vector grad_lam1(real_t x, real_t y, real_t z);
103 static Vector grad_lam2(real_t x, real_t y, real_t z);
104 static Vector grad_lam3(real_t x, real_t y, real_t z);
105 static Vector grad_lam4(real_t x, real_t y, real_t z);
106 static Vector grad_lam5(real_t x, real_t y, real_t z);
107
108 /// Two component vectors associated with edges touching the apex
110 { return Vector({lam1(x, y, z), lam5(x, y, z)}); }
112 { return Vector({lam2(x, y, z), lam5(x, y, z)}); }
114 { return Vector({lam3(x, y, z), lam5(x, y, z)}); }
116 { return Vector({lam4(x, y, z), lam5(x, y, z)}); }
117
118 /// Gradients of the above two component vectors
119 static DenseMatrix grad_lam15(real_t x, real_t y, real_t z);
120 static DenseMatrix grad_lam25(real_t x, real_t y, real_t z);
121 static DenseMatrix grad_lam35(real_t x, real_t y, real_t z);
122 static DenseMatrix grad_lam45(real_t x, real_t y, real_t z);
123
124 /// Computes $\lambda_i \nabla \lambda_5 - \lambda_5 \nabla \lambda_i$
129
130 /// Three component vectors associated with triangular faces
132 { return Vector({lam1(x, y, z), lam2(x, y, z), lam5(x, y, z)}); }
134 { return Vector({lam2(x, y, z), lam3(x, y, z), lam5(x, y, z)}); }
136 { return Vector({lam3(x, y, z), lam4(x, y, z), lam5(x, y, z)}); }
138 { return Vector({lam4(x, y, z), lam3(x, y, z), lam5(x, y, z)}); }
140 { return Vector({lam4(x, y, z), lam1(x, y, z), lam5(x, y, z)}); }
142 { return Vector({lam1(x, y, z), lam4(x, y, z), lam5(x, y, z)}); }
143
144 /// Vector functions related to the normals to the triangular faces
145 ///
146 /// Computes
147 /// $
148 /// \lambda_i \nabla\lambda_j \times \nabla \lambda_5
149 /// + \lambda_j \nabla\lambda_5 \times \nabla \lambda_i
150 /// + \lambda_5 \nabla\lambda_i \times \nabla \lambda_j
151 /// $
158
159 /// Divergences of the above "normal" vector functions divided by 3
166
167 static real_t mu0(real_t z)
168 { return 1.0 - z; }
169 static real_t mu1(real_t z)
170 { return z; }
171
173 { return Vector({0.0, 0.0, -1.0}); }
175 { return Vector({0.0, 0.0, 1.0}); }
176
178 { return Vector({mu0(z), mu1(z)}); }
179
180 static DenseMatrix grad_mu01(real_t z);
181
182 static real_t mu0(real_t z, const Vector &xy, unsigned int ab)
183 { return 1.0 - xy[ab-1] / (1.0 - z); }
184 static real_t mu1(real_t z, const Vector &xy, unsigned int ab)
185 { return xy[ab-1] / (1.0 - z); }
186
187 static Vector grad_mu0(real_t z, const Vector xy, unsigned int ab);
188 static Vector grad_mu1(real_t z, const Vector xy, unsigned int ab);
189
190 static Vector mu01(real_t z, Vector xy, unsigned int ab)
191 { return Vector({mu0(z, xy, ab), mu1(z, xy, ab)}); }
192
193 static DenseMatrix grad_mu01(real_t z, Vector xy, unsigned int ab);
194 static Vector mu01_grad_mu01(real_t z, Vector xy, unsigned int ab);
195
196 static real_t nu0(real_t z, Vector xy, unsigned int ab)
197 { return 1.0 - xy[ab-1] - z; }
198 static real_t nu1(real_t z, Vector xy, unsigned int ab) { return xy[ab-1]; }
199 static real_t nu2(real_t z, Vector xy, unsigned int ab) { return z; }
200
201 static Vector grad_nu0(real_t z, const Vector xy, unsigned int ab);
202 static Vector grad_nu1(real_t z, const Vector xy, unsigned int ab);
203 static Vector grad_nu2(real_t z, const Vector xy, unsigned int ab);
204
205 static Vector nu01(real_t z, Vector xy, unsigned int ab)
206 { return Vector({nu0(z, xy, ab), nu1(z, xy, ab)}); }
207 static Vector nu12(real_t z, Vector xy, unsigned int ab)
208 { return Vector({nu1(z, xy, ab), nu2(z, xy, ab)}); }
209 static Vector nu012(real_t z, Vector xy, unsigned int ab)
210 { return Vector({nu0(z, xy, ab), nu1(z, xy, ab), nu2(z, xy, ab)}); }
211 static Vector nu120(real_t z, Vector xy, unsigned int ab)
212 { return Vector({nu1(z, xy, ab), nu2(z, xy, ab), nu0(z, xy, ab)}); }
213
214 static DenseMatrix grad_nu01(real_t z, Vector xy, unsigned int ab);
215 static DenseMatrix grad_nu012(real_t z, Vector xy, unsigned int ab);
216 static DenseMatrix grad_nu120(real_t z, Vector xy, unsigned int ab);
217
218 static Vector nu01_grad_nu01(real_t z, Vector xy, unsigned int ab);
219 static Vector nu12_grad_nu12(real_t z, Vector xy, unsigned int ab);
220 static Vector nu012_grad_nu012(real_t z, Vector xy, unsigned int ab);
221
222 /// Shifted and Scaled Legendre Polynomials
223 /** Implements a scaled and shifted set of Legendre polynomials
224
225 $P_i(x;t) = P_i(x / t) * t^i$
226
227 where @a t >= 0.0, @a x $\in [0,t]$, and $P_i$ is the shifted Legendre
228 polynomial defined on $[0,1]$ rather than the usual $[-1,1]$. The
229 entries stored in @a u correspond to the values of
230 $P_0$, $P_1$, ... $P_p$.
231
232 @a u must be at least @a p + 1 in length
233 */
234 static void CalcScaledLegendre(int p, real_t x, real_t t,
235 real_t *u);
236 static void CalcScaledLegendre(int p, real_t x, real_t t,
237 real_t *u, real_t *dudx, real_t *dudt);
238
239 static void CalcScaledLegendre(int p, real_t x, real_t t,
240 Vector &u);
241 static void CalcScaledLegendre(int p, real_t x, real_t t,
242 Vector &u, Vector &dudx, Vector &dudt);
243
244 /// Integrated Legendre Polynomials
245 /** These are the integrals of the shifted and scaled Legendre polynomials
246 provided above and defined as:
247
248 $L_i(x;t) = \int_0^x P_{i-1}(y;t)dy\mbox{ for }i>=1$
249
250 These polynomials are computed as:
251
252 $L_0(x;t) = 0$, $L_1(x;t) = x$,
253
254 $2(2i-1)L_i(x;t) = P_i(x;t) - t^2 P_{i-2}(x;t)\mbox{ for }i>=2$
255
256 @a u must be at least @a p + 1 in length
257 */
258 static void CalcIntegratedLegendre(int p, real_t x,
259 real_t t, real_t *u);
260 static void CalcIntegratedLegendre(int p, real_t x,
261 real_t t, real_t *u,
262 real_t *dudx, real_t *dudt);
263
264 static void CalcIntegratedLegendre(int p, real_t x,
265 real_t t, Vector &u);
266 static void CalcIntegratedLegendre(int p, real_t x,
267 real_t t, Vector &u,
268 Vector &dudx, Vector &dudt);
269
270 /** @a u must be at least @a p + 1 in length */
271 static void CalcHomogenizedScaLegendre(int p, real_t s0, real_t s1,
272 real_t *u);
273 static void CalcHomogenizedScaLegendre(int p,
274 real_t s0, real_t s1,
275 real_t *u,
276 real_t *duds0, real_t *duds1);
277 static void CalcHomogenizedScaLegendre(int p, real_t s0, real_t s1,
278 Vector &u);
279 static void CalcHomogenizedScaLegendre(int p,
280 real_t s0, real_t s1,
281 Vector &u,
282 Vector &duds0, Vector &duds1);
283
284 /** @a u must be at least @a p + 1 in length */
285 static void CalcHomogenizedIntLegendre(int p,
286 real_t t0, real_t t1,
287 real_t *u);
288 static void CalcHomogenizedIntLegendre(int p,
289 real_t t0, real_t t1,
290 real_t *u,
291 real_t *dudt0, real_t *dudt1);
292 static void CalcHomogenizedIntLegendre(int p,
293 real_t t0, real_t t1,
294 Vector &u);
295 static void CalcHomogenizedIntLegendre(int p,
296 real_t t0, real_t t1,
297 Vector &u,
298 Vector &dudt0, Vector &dudt1);
299
300 /// Shifted and Scaled Jacobi Polynomials
301 /** Implements a scaled and shifted set of Jacobi polynomials
302
303 $P^\alpha_i(x / t) * t^i$
304
305 where @a alpha $= \alpha >-1$, @a t $>= 0.0$, @a x $\in [0,t]$, and
306 $P^\alpha_i$ is the shifted Jacobi polynomial defined on $[0,1]$ rather
307 than the usual $[-1,1]$. The entries stored in @a u correspond to the
308 values of $P^\alpha_0$, $P^\alpha_1$, ... $P^\alpha_p$.
309
310 @note Jacobi polynomials typically posses two parameters,
311 $P^{\alpha, \beta}_i$, but we only consider the special case where
312 $\beta=0$.
313
314 @a u must be at least @a p + 1 in length
315 */
316 static void CalcScaledJacobi(int p, real_t alpha,
317 real_t x, real_t t,
318 real_t *u);
319 static void CalcScaledJacobi(int p, real_t alpha,
320 real_t x, real_t t,
321 real_t *u, real_t *dudx, real_t *dudt);
322
323 static void CalcScaledJacobi(int p, real_t alpha,
324 real_t x, real_t t,
325 Vector &u);
326 static void CalcScaledJacobi(int p, real_t alpha,
327 real_t x, real_t t,
328 Vector &u, Vector &dudx, Vector &dudt);
329
330 /// Integrated Jacobi Polynomials
331 /** These are the integrals of the shifted and scaled Jacobi polynomials
332 provided above and defined as:
333
334 $L^\alpha_i(x;t) = \int_0^x P^\alpha_{i-1}(y;t)dy\mbox{ for }i>=1$
335
336 These polynomials are computed as:
337
338 $L^\alpha_0(x;t) = 0$, $L^\alpha_1(x;t) = x$,
339
340 $L^\alpha_i(x;t) = a_i P^\alpha_i(x;t) + b_i t P^\alpha_{i-1}(x;t)
341 - c_i t^2 P^\alpha_{i-2}(x;t)\mbox{ for }i>=2$
342
343 With
344
345 $a_i = (i + \alpha) / (2i + \alpha - 1)(2i + \alpha)$
346
347 $b_i = \alpha / (2i + \alpha - 2)(2i + \alpha)$
348
349 $c_i = (i - 1) / (2i + \alpha - 2)(2i + \alpha - 1)$
350
351 @a u must be at least @a p + 1 in length
352 */
353 static void CalcIntegratedJacobi(int p, real_t alpha,
354 real_t x, real_t t,
355 real_t *u);
356 static void CalcIntegratedJacobi(int p, real_t alpha,
357 real_t x, real_t t,
358 real_t *u, real_t *dudx, real_t *dudt);
359
361 real_t x, real_t t,
362 Vector &u)
363 { CalcIntegratedJacobi(p, alpha, x, t, u.GetData()); }
365 real_t x, real_t t,
366 Vector &u, Vector &dudx, Vector &dudt)
367 {
368 CalcIntegratedJacobi(p, alpha, x, t, u.GetData(),
369 dudx.GetData(), dudt.GetData());
370 }
371
372 /** @a u must be at least @a p + 1 in length */
374 real_t t0, real_t t1,
375 real_t *u)
376 { CalcScaledJacobi(p, alpha, t1, t0 + t1, u); }
377 static void CalcHomogenizedScaJacobi(int p, real_t alpha,
378 real_t t0, real_t t1,
379 real_t *u,
380 real_t *dudt0, real_t *dudt1);
381 static void CalcHomogenizedScaJacobi(int p, real_t alpha,
382 real_t t0, real_t t1,
383 Vector &u);
384 static void CalcHomogenizedScaJacobi(int p, real_t alpha,
385 real_t t0, real_t t1,
386 Vector &u,
387 Vector &dudt0, Vector &dudt1);
388
389 /** @a u must be at least @a p + 1 in length */
391 real_t t0, real_t t1,
392 real_t *u)
393 { CalcIntegratedJacobi(p, alpha, t1, t0 + t1, u); }
394 static void CalcHomogenizedIntJacobi(int p, real_t alpha,
395 real_t t0, real_t t1,
396 real_t *u,
397 real_t *dudt0, real_t *dudt1);
398 static void CalcHomogenizedIntJacobi(int p, real_t alpha,
399 real_t t0, real_t t1,
400 Vector &u);
401 static void CalcHomogenizedIntJacobi(int p, real_t alpha,
402 real_t t0, real_t t1,
403 Vector &u,
404 Vector &dudt0, Vector &dudt1);
405
406 /** @a u must be at least @a p + 1 in length */
407 static void phi_E(int p, real_t s0, real_t s1, real_t *u);
408 static void phi_E(int p, real_t s0, real_t s1, real_t *u,
409 real_t *duds0, real_t *duds1);
410 static void phi_E(int p, Vector s, Vector &u);
411 static void phi_E(int p, Vector s, Vector &u, DenseMatrix &duds);
412
413 /** @a grad_s must be 2x3 */
414 void phi_E(int p, Vector s, const DenseMatrix &grad_s,
415 Vector &u, DenseMatrix &grad_u) const;
416
417 /** @a u must be at least (p+1)x(p+1) in size */
418 void phi_Q(int p, Vector s, Vector t, DenseMatrix &u) const;
419 void phi_Q(int p, Vector s, const DenseMatrix &grad_s,
420 Vector t, const DenseMatrix &grad_t,
421 DenseMatrix &u, DenseTensor &grad_u) const;
422 void phi_T(int p, Vector nu, DenseMatrix &u) const;
423 void phi_T(int p, Vector nu, const DenseMatrix &grad_nu,
424 DenseMatrix &u, DenseTensor &grad_u) const;
425
426 /** This is a vector-valued function associated with an edge of a pyramid
427
428 The vector @a s contains two coordinate values and @a ds is related to the
429 gradient of these coordinates with respect to the reference coordinates i.e.
430 sds = s0 grad s1 - s1 grad s0
431 */
432 void E_E(int p, Vector s, Vector sds, DenseMatrix &u) const;
433 void E_E(int p, Vector s, const DenseMatrix &grad_s, DenseMatrix &u,
434 DenseMatrix &curl_u) const;
435
436 void E_Q(int p, Vector s, Vector ds, Vector t,
437 DenseTensor &u) const;
438 void E_Q(int p, Vector s, const DenseMatrix &grad_s,
439 Vector t, const DenseMatrix &grad_t,
440 DenseTensor &u, DenseTensor &curl_u) const;
441
442 void E_T(int p, Vector s, Vector sds, DenseTensor &u) const;
443 void E_T(int p, Vector s, const DenseMatrix &grad_s,
444 DenseTensor &u, DenseTensor &curl_u) const;
445
446 /** This is a vector-valued function associated with the quadrilateral face
447 of a pyramid
448
449 The vectors @a s and @a t contain pairs of coordinate values and @a ds and
450 @a dt are related to derivatives of these coordinates:
451
452 ds = s0 grad s1 - s1 grad s0
453
454 dt = t0 grad t1 - t1 grad t0
455 */
456 void V_Q(int p, Vector s, Vector ds, Vector t, Vector dt,
457 DenseTensor &u) const;
458
459 /** This is a vector-valued function associated with the triangular faces of
460 a pyramid
461
462 The vector @a s contains three coordinate values and @a sdsxds is related to
463 derivatives of these coordinates with respect to the reference coordinates:
464
465 sdsxds = s0 grad s1 x grad s2 + s1 grad s2 x grad s0 +
466 s2 grad s0 x grad s1
467 */
468 void V_T(int p, Vector s, Vector sdsxds, DenseTensor &u) const;
469
470 /** This computes V_T as above and its divergence
471
472 The vector @a s contains three coordinate values and @a sdsxds is related to
473 derivatives of these coordinates with respect to the reference coordinates:
474
475 sdsxds = s0 grad s1 x grad s2 + s1 grad s2 x grad s0 +
476 s2 grad s0 x grad s1
477
478 The scalar @a dsdsxds is the divergence of sdsxds:
479
480 dsdsxds = grad s0 dot (grad s1 x grad s2)
481 */
482 void V_T(int p, Vector s, Vector sdsxds, real_t dsdsxds,
483 DenseTensor &u, DenseMatrix &du) const;
484
485 void VT_T(int p, Vector s, Vector sds, Vector sdsxds,
486 real_t mu, Vector grad_mu, DenseTensor &u) const;
487 void VT_T(int p, Vector s, Vector sds, Vector sdsxds,
488 Vector grad_s2, real_t mu, Vector grad_mu,
489 DenseTensor &u, DenseMatrix &du) const;
490
491 /** This implements $V^\unlhd_{ij}$ from the Fuentes paper
492
493 @a u must be at least (p+1)x(p+1)x3
494 */
495 void V_L(int p, Vector sx, const DenseMatrix &grad_sx,
496 Vector sy, const DenseMatrix &grad_sy,
497 real_t t, Vector grad_t, DenseTensor &u) const;
498
499 /** This implements $V^\unrhd_i$ from the Fuentes paper
500
501 @a u must be at least (p+1)x3 */
502 void V_R(int p, Vector s, const DenseMatrix &grad_s,
503 real_t mu, Vector grad_mu,
504 real_t t, Vector grad_t, DenseMatrix &u) const;
505};
506
507} // namespace mfem
508
509#endif
510
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Rank 3 tensor (array of matrices)
static DenseMatrix grad_mu01(real_t z)
static Vector lam45_grad_lam45(real_t x, real_t y, real_t z)
void phi_Q(int p, Vector s, Vector t, DenseMatrix &u) const
static void CalcHomogenizedIntJacobi(int p, real_t alpha, real_t t0, real_t t1, real_t *u)
static Vector nu12(real_t z, Vector xy, unsigned int ab)
void V_R(int p, Vector s, const DenseMatrix &grad_s, real_t mu, Vector grad_mu, real_t t, Vector grad_t, DenseMatrix &u) const
static Vector lam125_grad_lam125(real_t x, real_t y, real_t z)
static constexpr real_t apex_tol
static Vector lam145_grad_lam145(real_t x, real_t y, real_t z)
static void CalcIntegratedJacobi(int p, real_t alpha, real_t x, real_t t, Vector &u)
static Vector lam45(real_t x, real_t y, real_t z)
static void phi_E(int p, real_t s0, real_t s1, real_t *u)
static Vector nu012(real_t z, Vector xy, unsigned int ab)
static DenseMatrix grad_nu120(real_t z, Vector xy, unsigned int ab)
static real_t lam4(real_t x, real_t y, real_t z)
static DenseMatrix grad_lam25(real_t x, real_t y, real_t z)
void E_E(int p, Vector s, Vector sds, DenseMatrix &u) const
static void CalcIntegratedLegendre(int p, real_t x, real_t t, real_t *u)
Integrated Legendre Polynomials.
static Vector nu01_grad_nu01(real_t z, Vector xy, unsigned int ab)
static Vector lam35(real_t x, real_t y, real_t z)
static Vector grad_lam5(real_t x, real_t y, real_t z)
static real_t lam3(real_t x, real_t y, real_t z)
static Vector lam235(real_t x, real_t y, real_t z)
static real_t div_lam345_grad_lam345(real_t x, real_t y, real_t z)
static Vector lam25(real_t x, real_t y, real_t z)
static Vector nu01(real_t z, Vector xy, unsigned int ab)
void VT_T(int p, Vector s, Vector sds, Vector sdsxds, real_t mu, Vector grad_mu, DenseTensor &u) const
static Vector lam25_grad_lam25(real_t x, real_t y, real_t z)
static Vector lam415_grad_lam415(real_t x, real_t y, real_t z)
static Vector lam235_grad_lam235(real_t x, real_t y, real_t z)
static Vector mu01(real_t z)
static DenseMatrix grad_nu012(real_t z, Vector xy, unsigned int ab)
static void CalcIntegratedJacobi(int p, real_t alpha, real_t x, real_t t, real_t *u)
Integrated Jacobi Polynomials.
FuentesPyramid()=default
static real_t div_lam145_grad_lam145(real_t x, real_t y, real_t z)
void V_L(int p, Vector sx, const DenseMatrix &grad_sx, Vector sy, const DenseMatrix &grad_sy, real_t t, Vector grad_t, DenseTensor &u) const
static Vector lam125(real_t x, real_t y, real_t z)
Three component vectors associated with triangular faces.
static Vector grad_lam4(real_t x, real_t y, real_t z)
static Vector grad_nu2(real_t z, const Vector xy, unsigned int ab)
static Vector grad_mu1(real_t z)
void V_Q(int p, Vector s, Vector ds, Vector t, Vector dt, DenseTensor &u) const
static void CalcScaledLegendre(int p, real_t x, real_t t, real_t *u)
Shifted and Scaled Legendre Polynomials.
void E_T(int p, Vector s, Vector sds, DenseTensor &u) const
static real_t mu0(real_t z)
static Vector mu01(real_t z, Vector xy, unsigned int ab)
static Vector nu120(real_t z, Vector xy, unsigned int ab)
static Vector lam345(real_t x, real_t y, real_t z)
static real_t lam2(real_t x, real_t y, real_t z)
static void CalcHomogenizedScaJacobi(int p, real_t alpha, real_t t0, real_t t1, real_t *u)
static Vector lam15_grad_lam15(real_t x, real_t y, real_t z)
Computes .
static real_t mu0(real_t z, const Vector &xy, unsigned int ab)
static Vector lam35_grad_lam35(real_t x, real_t y, real_t z)
static real_t div_lam415_grad_lam415(real_t x, real_t y, real_t z)
static real_t lam1(real_t x, real_t y, real_t z)
Pyramid "Affine" Coordinates.
void V_T(int p, Vector s, Vector sdsxds, DenseTensor &u) const
static Vector lam435_grad_lam435(real_t x, real_t y, real_t z)
static real_t div_lam235_grad_lam235(real_t x, real_t y, real_t z)
static Vector grad_mu0(real_t z)
static constexpr real_t one
static DenseMatrix grad_lam35(real_t x, real_t y, real_t z)
static real_t lam5(real_t x, real_t y, real_t z)
static real_t mu1(real_t z, const Vector &xy, unsigned int ab)
static Vector lam415(real_t x, real_t y, real_t z)
static real_t nu0(real_t z, Vector xy, unsigned int ab)
static DenseMatrix grad_lam15(real_t x, real_t y, real_t z)
Gradients of the above two component vectors.
static Vector grad_nu1(real_t z, const Vector xy, unsigned int ab)
void E_Q(int p, Vector s, Vector ds, Vector t, DenseTensor &u) const
static Vector grad_lam1(real_t x, real_t y, real_t z)
Gradients of the "Affine" Coordinates.
static Vector lam15(real_t x, real_t y, real_t z)
Two component vectors associated with edges touching the apex.
static Vector mu01_grad_mu01(real_t z, Vector xy, unsigned int ab)
static Vector lam145(real_t x, real_t y, real_t z)
static real_t mu1(real_t z)
static real_t nu2(real_t z, Vector xy, unsigned int ab)
static Vector nu012_grad_nu012(real_t z, Vector xy, unsigned int ab)
static Vector grad_lam2(real_t x, real_t y, real_t z)
static Vector lam345_grad_lam345(real_t x, real_t y, real_t z)
static void CalcHomogenizedIntLegendre(int p, real_t t0, real_t t1, real_t *u)
static real_t div_lam125_grad_lam125(real_t x, real_t y, real_t z)
Divergences of the above "normal" vector functions divided by 3.
static DenseMatrix grad_lam45(real_t x, real_t y, real_t z)
static DenseMatrix grad_nu01(real_t z, Vector xy, unsigned int ab)
static real_t nu1(real_t z, Vector xy, unsigned int ab)
static Vector lam435(real_t x, real_t y, real_t z)
static void CalcIntegratedJacobi(int p, real_t alpha, real_t x, real_t t, Vector &u, Vector &dudx, Vector &dudt)
static bool CheckZ(real_t z)
static void CalcScaledJacobi(int p, real_t alpha, real_t x, real_t t, real_t *u)
Shifted and Scaled Jacobi Polynomials.
static void CalcHomogenizedScaLegendre(int p, real_t s0, real_t s1, real_t *u)
static Vector grad_nu0(real_t z, const Vector xy, unsigned int ab)
void phi_T(int p, Vector nu, DenseMatrix &u) const
static Vector nu12_grad_nu12(real_t z, Vector xy, unsigned int ab)
static real_t div_lam435_grad_lam435(real_t x, real_t y, real_t z)
static Vector grad_lam3(real_t x, real_t y, real_t z)
static constexpr real_t zero
Vector data type.
Definition vector.hpp:82
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:235
const real_t alpha
Definition ex15.cpp:369
real_t mu
Definition ex25.cpp:140
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
float real_t
Definition config.hpp:43
real_t p(const Vector &x, real_t t)