MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::FuentesPyramid Class Reference

#include <fe_pyramid.hpp>

Inheritance diagram for mfem::FuentesPyramid:
[legend]

Public Member Functions

 FuentesPyramid ()=default
 
void phi_E (int p, Vector s, const DenseMatrix &grad_s, Vector &u, DenseMatrix &grad_u) const
 
void phi_Q (int p, Vector s, Vector t, DenseMatrix &u) const
 
void phi_Q (int p, Vector s, const DenseMatrix &grad_s, Vector t, const DenseMatrix &grad_t, DenseMatrix &u, DenseTensor &grad_u) const
 
void phi_T (int p, Vector nu, DenseMatrix &u) const
 
void phi_T (int p, Vector nu, const DenseMatrix &grad_nu, DenseMatrix &u, DenseTensor &grad_u) const
 
void E_E (int p, Vector s, Vector sds, DenseMatrix &u) const
 
void E_E (int p, Vector s, const DenseMatrix &grad_s, DenseMatrix &u, DenseMatrix &curl_u) const
 
void E_Q (int p, Vector s, Vector ds, Vector t, DenseTensor &u) const
 
void E_Q (int p, Vector s, const DenseMatrix &grad_s, Vector t, const DenseMatrix &grad_t, DenseTensor &u, DenseTensor &curl_u) const
 
void E_T (int p, Vector s, Vector sds, DenseTensor &u) const
 
void E_T (int p, Vector s, const DenseMatrix &grad_s, DenseTensor &u, DenseTensor &curl_u) const
 
void V_Q (int p, Vector s, Vector ds, Vector t, Vector dt, DenseTensor &u) const
 
void V_T (int p, Vector s, Vector sdsxds, DenseTensor &u) const
 
void V_T (int p, Vector s, Vector sdsxds, real_t dsdsxds, DenseTensor &u, DenseMatrix &du) const
 
void VT_T (int p, Vector s, Vector sds, Vector sdsxds, real_t mu, Vector grad_mu, DenseTensor &u) const
 
void VT_T (int p, Vector s, Vector sds, Vector sdsxds, Vector grad_s2, real_t mu, Vector grad_mu, DenseTensor &u, DenseMatrix &du) const
 
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
 
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 Public Member Functions

static bool CheckZ (real_t z)
 
static real_t lam1 (real_t x, real_t y, real_t z)
 Pyramid "Affine" Coordinates.
 
static real_t lam2 (real_t x, real_t y, real_t z)
 
static real_t lam3 (real_t x, real_t y, real_t z)
 
static real_t lam4 (real_t x, real_t y, real_t z)
 
static real_t lam5 (real_t x, real_t y, real_t z)
 
static Vector grad_lam1 (real_t x, real_t y, real_t z)
 Gradients of the "Affine" Coordinates.
 
static Vector grad_lam2 (real_t x, real_t y, real_t z)
 
static Vector grad_lam3 (real_t x, real_t y, real_t z)
 
static Vector grad_lam4 (real_t x, real_t y, real_t z)
 
static Vector grad_lam5 (real_t x, real_t y, real_t z)
 
static Vector lam15 (real_t x, real_t y, real_t z)
 Two component vectors associated with edges touching the apex.
 
static Vector lam25 (real_t x, real_t y, real_t z)
 
static Vector lam35 (real_t x, real_t y, real_t z)
 
static Vector lam45 (real_t x, real_t y, real_t z)
 
static DenseMatrix grad_lam15 (real_t x, real_t y, real_t z)
 Gradients of the above two component vectors.
 
static DenseMatrix grad_lam25 (real_t x, real_t y, real_t z)
 
static DenseMatrix grad_lam35 (real_t x, real_t y, real_t z)
 
static DenseMatrix grad_lam45 (real_t x, real_t y, real_t z)
 
static Vector lam15_grad_lam15 (real_t x, real_t y, real_t z)
 Computes λiλ5λ5λi.
 
static Vector lam25_grad_lam25 (real_t x, real_t y, real_t z)
 
static Vector lam35_grad_lam35 (real_t x, real_t y, real_t z)
 
static Vector lam45_grad_lam45 (real_t x, real_t y, real_t z)
 
static Vector lam125 (real_t x, real_t y, real_t z)
 Three component vectors associated with triangular faces.
 
static Vector lam235 (real_t x, real_t y, real_t z)
 
static Vector lam345 (real_t x, real_t y, real_t z)
 
static Vector lam435 (real_t x, real_t y, real_t z)
 
static Vector lam415 (real_t x, real_t y, real_t z)
 
static Vector lam145 (real_t x, real_t y, real_t z)
 
static Vector lam125_grad_lam125 (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 lam345_grad_lam345 (real_t x, real_t y, real_t z)
 
static Vector lam435_grad_lam435 (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 lam145_grad_lam145 (real_t x, real_t y, real_t z)
 
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 real_t div_lam235_grad_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 real_t div_lam435_grad_lam435 (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 div_lam145_grad_lam145 (real_t x, real_t y, real_t z)
 
static real_t mu0 (real_t z)
 
static real_t mu1 (real_t z)
 
static Vector grad_mu0 (real_t z)
 
static Vector grad_mu1 (real_t z)
 
static Vector mu01 (real_t z)
 
static DenseMatrix grad_mu01 (real_t z)
 
static real_t mu0 (real_t z, const Vector &xy, unsigned int ab)
 
static real_t mu1 (real_t z, const Vector &xy, unsigned int ab)
 
static Vector grad_mu0 (real_t z, const Vector xy, unsigned int ab)
 
static Vector grad_mu1 (real_t z, const Vector xy, unsigned int ab)
 
static Vector mu01 (real_t z, Vector xy, unsigned int ab)
 
static DenseMatrix grad_mu01 (real_t z, Vector xy, unsigned int ab)
 
static Vector mu01_grad_mu01 (real_t z, Vector xy, unsigned int ab)
 
static real_t nu0 (real_t z, Vector xy, unsigned int ab)
 
static real_t nu1 (real_t z, Vector xy, unsigned int ab)
 
static real_t nu2 (real_t z, Vector xy, unsigned int ab)
 
static Vector grad_nu0 (real_t z, const Vector xy, unsigned int ab)
 
static Vector grad_nu1 (real_t z, const Vector xy, unsigned int ab)
 
static Vector grad_nu2 (real_t z, const Vector xy, unsigned int ab)
 
static Vector nu01 (real_t z, Vector xy, unsigned int ab)
 
static Vector nu12 (real_t z, Vector xy, unsigned int ab)
 
static Vector nu012 (real_t z, Vector xy, unsigned int ab)
 
static Vector nu120 (real_t z, Vector xy, unsigned int ab)
 
static DenseMatrix grad_nu01 (real_t z, Vector xy, unsigned int ab)
 
static DenseMatrix grad_nu012 (real_t z, Vector xy, unsigned int ab)
 
static DenseMatrix grad_nu120 (real_t z, Vector xy, unsigned int ab)
 
static Vector nu01_grad_nu01 (real_t z, Vector xy, unsigned int ab)
 
static Vector nu12_grad_nu12 (real_t z, Vector xy, unsigned int ab)
 
static Vector nu012_grad_nu012 (real_t z, Vector xy, unsigned int ab)
 
static void CalcScaledLegendre (int p, real_t x, real_t t, real_t *u)
 Shifted and Scaled Legendre Polynomials.
 
static void CalcScaledLegendre (int p, real_t x, real_t t, real_t *u, real_t *dudx, real_t *dudt)
 
static void CalcScaledLegendre (int p, real_t x, real_t t, Vector &u)
 
static void CalcScaledLegendre (int p, real_t x, real_t t, Vector &u, Vector &dudx, Vector &dudt)
 
static void CalcIntegratedLegendre (int p, real_t x, real_t t, real_t *u)
 Integrated Legendre Polynomials.
 
static void CalcIntegratedLegendre (int p, real_t x, real_t t, real_t *u, real_t *dudx, real_t *dudt)
 
static void CalcIntegratedLegendre (int p, real_t x, real_t t, Vector &u)
 
static void CalcIntegratedLegendre (int p, real_t x, real_t t, Vector &u, Vector &dudx, Vector &dudt)
 
static void CalcHomogenizedScaLegendre (int p, real_t s0, real_t s1, real_t *u)
 
static void CalcHomogenizedScaLegendre (int p, real_t s0, real_t s1, real_t *u, real_t *duds0, real_t *duds1)
 
static void CalcHomogenizedScaLegendre (int p, real_t s0, real_t s1, Vector &u)
 
static void CalcHomogenizedScaLegendre (int p, real_t s0, real_t s1, Vector &u, Vector &duds0, Vector &duds1)
 
static void CalcHomogenizedIntLegendre (int p, real_t t0, real_t t1, real_t *u)
 
static void CalcHomogenizedIntLegendre (int p, real_t t0, real_t t1, real_t *u, real_t *dudt0, real_t *dudt1)
 
static void CalcHomogenizedIntLegendre (int p, real_t t0, real_t t1, Vector &u)
 
static void CalcHomogenizedIntLegendre (int p, real_t t0, real_t t1, Vector &u, Vector &dudt0, Vector &dudt1)
 
static void CalcScaledJacobi (int p, real_t alpha, real_t x, real_t t, real_t *u)
 Shifted and Scaled Jacobi Polynomials.
 
static void CalcScaledJacobi (int p, real_t alpha, real_t x, real_t t, real_t *u, real_t *dudx, real_t *dudt)
 
static void CalcScaledJacobi (int p, real_t alpha, real_t x, real_t t, Vector &u)
 
static void CalcScaledJacobi (int p, real_t alpha, real_t x, real_t t, Vector &u, Vector &dudx, Vector &dudt)
 
static void CalcIntegratedJacobi (int p, real_t alpha, real_t x, real_t t, real_t *u)
 Integrated Jacobi Polynomials.
 
static void CalcIntegratedJacobi (int p, real_t alpha, real_t x, real_t t, real_t *u, real_t *dudx, real_t *dudt)
 
static void CalcIntegratedJacobi (int p, real_t alpha, real_t x, real_t t, Vector &u)
 
static void CalcIntegratedJacobi (int p, real_t alpha, real_t x, real_t t, Vector &u, Vector &dudx, Vector &dudt)
 
static void CalcHomogenizedScaJacobi (int p, real_t alpha, real_t t0, real_t t1, real_t *u)
 
static void CalcHomogenizedScaJacobi (int p, real_t alpha, real_t t0, real_t t1, real_t *u, real_t *dudt0, real_t *dudt1)
 
static void CalcHomogenizedScaJacobi (int p, real_t alpha, real_t t0, real_t t1, Vector &u)
 
static void CalcHomogenizedScaJacobi (int p, real_t alpha, real_t t0, real_t t1, Vector &u, Vector &dudt0, Vector &dudt1)
 
static void CalcHomogenizedIntJacobi (int p, real_t alpha, real_t t0, real_t t1, real_t *u)
 
static void CalcHomogenizedIntJacobi (int p, real_t alpha, real_t t0, real_t t1, real_t *u, real_t *dudt0, real_t *dudt1)
 
static void CalcHomogenizedIntJacobi (int p, real_t alpha, real_t t0, real_t t1, Vector &u)
 
static void CalcHomogenizedIntJacobi (int p, real_t alpha, real_t t0, real_t t1, Vector &u, Vector &dudt0, Vector &dudt1)
 
static void phi_E (int p, real_t s0, real_t s1, real_t *u)
 
static void phi_E (int p, real_t s0, real_t s1, real_t *u, real_t *duds0, real_t *duds1)
 
static void phi_E (int p, Vector s, Vector &u)
 
static void phi_E (int p, Vector s, Vector &u, DenseMatrix &duds)
 

Static Protected Attributes

static constexpr real_t one = 1.0
 
static constexpr real_t zero = 0.0
 
static constexpr real_t apex_tol = 1e-8
 

Detailed Description

Base class for arbitrary order basis functions on pyramid-shaped elements

This base class provides a common class to store temporary vectors, matrices, and tensors computed by various functions defined on pyramid-shaped elements.

The function names defined here are chosen to reflect, as closely as possible, those used in the paper "Orientation embedded high order shape functions for the exact sequence elements of all shapes" by Federico Fuentes, Brendan Keith, Leszek Demkowicz, and Sriram Nagaraj, see https://doi.org/10.1016/j.camwa.2015.04.027.

Note
Many of the functions below, e.g. lam1, lam2, etc. and related functions, are singular or multi-valued at the apex of the pyramid. The values returned near the apex are computed in the limit z->1 using (x, y, z) = ((1-z)/2, (1-z)/2, z) i.e. along the line from the center of the base to the apex.

Definition at line 38 of file fe_pyramid.hpp.

Constructor & Destructor Documentation

◆ FuentesPyramid()

mfem::FuentesPyramid::FuentesPyramid ( )
default

Member Function Documentation

◆ CalcHomogenizedIntJacobi() [1/4]

static void mfem::FuentesPyramid::CalcHomogenizedIntJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
real_t * u )
inlinestatic

u must be at least p + 1 in length

Definition at line 390 of file fe_pyramid.hpp.

◆ CalcHomogenizedIntJacobi() [2/4]

void mfem::FuentesPyramid::CalcHomogenizedIntJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
real_t * u,
real_t * dudt0,
real_t * dudt1 )
static

Definition at line 698 of file fe_pyramid.cpp.

◆ CalcHomogenizedIntJacobi() [3/4]

void mfem::FuentesPyramid::CalcHomogenizedIntJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
Vector & u )
static

Definition at line 708 of file fe_pyramid.cpp.

◆ CalcHomogenizedIntJacobi() [4/4]

void mfem::FuentesPyramid::CalcHomogenizedIntJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
Vector & u,
Vector & dudt0,
Vector & dudt1 )
static

Definition at line 717 of file fe_pyramid.cpp.

◆ CalcHomogenizedIntLegendre() [1/4]

void mfem::FuentesPyramid::CalcHomogenizedIntLegendre ( int p,
real_t t0,
real_t t1,
real_t * u )
static

u must be at least p + 1 in length

Definition at line 485 of file fe_pyramid.cpp.

◆ CalcHomogenizedIntLegendre() [2/4]

void mfem::FuentesPyramid::CalcHomogenizedIntLegendre ( int p,
real_t t0,
real_t t1,
real_t * u,
real_t * dudt0,
real_t * dudt1 )
static

Definition at line 493 of file fe_pyramid.cpp.

◆ CalcHomogenizedIntLegendre() [3/4]

void mfem::FuentesPyramid::CalcHomogenizedIntLegendre ( int p,
real_t t0,
real_t t1,
Vector & u )
static

Definition at line 503 of file fe_pyramid.cpp.

◆ CalcHomogenizedIntLegendre() [4/4]

void mfem::FuentesPyramid::CalcHomogenizedIntLegendre ( int p,
real_t t0,
real_t t1,
Vector & u,
Vector & dudt0,
Vector & dudt1 )
static

Definition at line 512 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaJacobi() [1/4]

static void mfem::FuentesPyramid::CalcHomogenizedScaJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
real_t * u )
inlinestatic

u must be at least p + 1 in length

Definition at line 373 of file fe_pyramid.hpp.

◆ CalcHomogenizedScaJacobi() [2/4]

void mfem::FuentesPyramid::CalcHomogenizedScaJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
real_t * u,
real_t * dudt0,
real_t * dudt1 )
static

Definition at line 666 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaJacobi() [3/4]

void mfem::FuentesPyramid::CalcHomogenizedScaJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
Vector & u )
static

Definition at line 676 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaJacobi() [4/4]

void mfem::FuentesPyramid::CalcHomogenizedScaJacobi ( int p,
real_t alpha,
real_t t0,
real_t t1,
Vector & u,
Vector & dudt0,
Vector & dudt1 )
static

Definition at line 685 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaLegendre() [1/4]

void mfem::FuentesPyramid::CalcHomogenizedScaLegendre ( int p,
real_t s0,
real_t s1,
real_t * u )
static

u must be at least p + 1 in length

Definition at line 447 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaLegendre() [2/4]

void mfem::FuentesPyramid::CalcHomogenizedScaLegendre ( int p,
real_t s0,
real_t s1,
real_t * u,
real_t * duds0,
real_t * duds1 )
static

Definition at line 454 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaLegendre() [3/4]

void mfem::FuentesPyramid::CalcHomogenizedScaLegendre ( int p,
real_t s0,
real_t s1,
Vector & u )
static

Definition at line 464 of file fe_pyramid.cpp.

◆ CalcHomogenizedScaLegendre() [4/4]

void mfem::FuentesPyramid::CalcHomogenizedScaLegendre ( int p,
real_t s0,
real_t s1,
Vector & u,
Vector & duds0,
Vector & duds1 )
static

Definition at line 472 of file fe_pyramid.cpp.

◆ CalcIntegratedJacobi() [1/4]

void mfem::FuentesPyramid::CalcIntegratedJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
real_t * u )
static

Integrated Jacobi Polynomials.

These are the integrals of the shifted and scaled Jacobi polynomials provided above and defined as:

Liα(x;t)=0xPi1α(y;t)dy for i>=1

These polynomials are computed as:

L0α(x;t)=0, L1α(x;t)=x,

Liα(x;t)=aiPiα(x;t)+bitPi1α(x;t)cit2Pi2α(x;t) for i>=2

With

ai=(i+α)/(2i+α1)(2i+α)

bi=α/(2i+α2)(2i+α)

ci=(i1)/(2i+α2)(2i+α1)

u must be at least p + 1 in length

Definition at line 525 of file fe_pyramid.cpp.

◆ CalcIntegratedJacobi() [2/4]

void mfem::FuentesPyramid::CalcIntegratedJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
real_t * u,
real_t * dudx,
real_t * dudt )
static

Definition at line 559 of file fe_pyramid.cpp.

◆ CalcIntegratedJacobi() [3/4]

static void mfem::FuentesPyramid::CalcIntegratedJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
Vector & u )
inlinestatic

Definition at line 360 of file fe_pyramid.hpp.

◆ CalcIntegratedJacobi() [4/4]

static void mfem::FuentesPyramid::CalcIntegratedJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
Vector & u,
Vector & dudx,
Vector & dudt )
inlinestatic

Definition at line 364 of file fe_pyramid.hpp.

◆ CalcIntegratedLegendre() [1/4]

void mfem::FuentesPyramid::CalcIntegratedLegendre ( int p,
real_t x,
real_t t,
real_t * u )
static

Integrated Legendre Polynomials.

These are the integrals of the shifted and scaled Legendre polynomials provided above and defined as:

Li(x;t)=0xPi1(y;t)dy for i>=1

These polynomials are computed as:

L0(x;t)=0, L1(x;t)=x,

2(2i1)Li(x;t)=Pi(x;t)t2Pi2(x;t) for i>=2

u must be at least p + 1 in length

Definition at line 369 of file fe_pyramid.cpp.

◆ CalcIntegratedLegendre() [2/4]

void mfem::FuentesPyramid::CalcIntegratedLegendre ( int p,
real_t x,
real_t t,
real_t * u,
real_t * dudx,
real_t * dudt )
static

Definition at line 395 of file fe_pyramid.cpp.

◆ CalcIntegratedLegendre() [3/4]

void mfem::FuentesPyramid::CalcIntegratedLegendre ( int p,
real_t x,
real_t t,
Vector & u )
static

Definition at line 427 of file fe_pyramid.cpp.

◆ CalcIntegratedLegendre() [4/4]

void mfem::FuentesPyramid::CalcIntegratedLegendre ( int p,
real_t x,
real_t t,
Vector & u,
Vector & dudx,
Vector & dudt )
static

Definition at line 435 of file fe_pyramid.cpp.

◆ CalcScaledJacobi() [1/4]

void mfem::FuentesPyramid::CalcScaledJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
real_t * u )
static

Shifted and Scaled Jacobi Polynomials.

Implements a scaled and shifted set of Jacobi polynomials

Piα(x/t)ti

where alpha =α>1, t >=0.0, x [0,t], and Piα is the shifted Jacobi polynomial defined on [0,1] rather than the usual [1,1]. The entries stored in u correspond to the values of P0α, P1α, ... Ppα.

Note
Jacobi polynomials typically posses two parameters, Piα,β, but we only consider the special case where β=0.

u must be at least p + 1 in length

Definition at line 591 of file fe_pyramid.cpp.

◆ CalcScaledJacobi() [2/4]

void mfem::FuentesPyramid::CalcScaledJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
real_t * u,
real_t * dudx,
real_t * dudt )
static

Definition at line 613 of file fe_pyramid.cpp.

◆ CalcScaledJacobi() [3/4]

void mfem::FuentesPyramid::CalcScaledJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
Vector & u )
static

Definition at line 645 of file fe_pyramid.cpp.

◆ CalcScaledJacobi() [4/4]

void mfem::FuentesPyramid::CalcScaledJacobi ( int p,
real_t alpha,
real_t x,
real_t t,
Vector & u,
Vector & dudx,
Vector & dudt )
static

Definition at line 654 of file fe_pyramid.cpp.

◆ CalcScaledLegendre() [1/4]

void mfem::FuentesPyramid::CalcScaledLegendre ( int p,
real_t x,
real_t t,
real_t * u )
static

Shifted and Scaled Legendre Polynomials.

Implements a scaled and shifted set of Legendre polynomials

Pi(x;t)=Pi(x/t)ti

where t >= 0.0, x [0,t], and Pi is the shifted Legendre polynomial defined on [0,1] rather than the usual [1,1]. The entries stored in u correspond to the values of P0, P1, ... Pp.

u must be at least p + 1 in length

Definition at line 293 of file fe_pyramid.cpp.

◆ CalcScaledLegendre() [2/4]

void mfem::FuentesPyramid::CalcScaledLegendre ( int p,
real_t x,
real_t t,
real_t * u,
real_t * dudx,
real_t * dudt )
static

Definition at line 313 of file fe_pyramid.cpp.

◆ CalcScaledLegendre() [3/4]

void mfem::FuentesPyramid::CalcScaledLegendre ( int p,
real_t x,
real_t t,
Vector & u )
static

Definition at line 351 of file fe_pyramid.cpp.

◆ CalcScaledLegendre() [4/4]

void mfem::FuentesPyramid::CalcScaledLegendre ( int p,
real_t x,
real_t t,
Vector & u,
Vector & dudx,
Vector & dudt )
static

Definition at line 359 of file fe_pyramid.cpp.

◆ CheckZ()

static bool mfem::FuentesPyramid::CheckZ ( real_t z)
inlinestatic

Definition at line 87 of file fe_pyramid.hpp.

◆ div_lam125_grad_lam125()

real_t mfem::FuentesPyramid::div_lam125_grad_lam125 ( real_t x,
real_t y,
real_t z )
static

Divergences of the above "normal" vector functions divided by 3.

Definition at line 160 of file fe_pyramid.cpp.

◆ div_lam145_grad_lam145()

real_t mfem::FuentesPyramid::div_lam145_grad_lam145 ( real_t x,
real_t y,
real_t z )
static

Definition at line 175 of file fe_pyramid.cpp.

◆ div_lam235_grad_lam235()

real_t mfem::FuentesPyramid::div_lam235_grad_lam235 ( real_t x,
real_t y,
real_t z )
static

Definition at line 163 of file fe_pyramid.cpp.

◆ div_lam345_grad_lam345()

real_t mfem::FuentesPyramid::div_lam345_grad_lam345 ( real_t x,
real_t y,
real_t z )
static

Definition at line 166 of file fe_pyramid.cpp.

◆ div_lam415_grad_lam415()

real_t mfem::FuentesPyramid::div_lam415_grad_lam415 ( real_t x,
real_t y,
real_t z )
static

Definition at line 172 of file fe_pyramid.cpp.

◆ div_lam435_grad_lam435()

real_t mfem::FuentesPyramid::div_lam435_grad_lam435 ( real_t x,
real_t y,
real_t z )
static

Definition at line 169 of file fe_pyramid.cpp.

◆ E_E() [1/2]

void mfem::FuentesPyramid::E_E ( int p,
Vector s,
const DenseMatrix & grad_s,
DenseMatrix & u,
DenseMatrix & curl_u ) const

Definition at line 983 of file fe_pyramid.cpp.

◆ E_E() [2/2]

void mfem::FuentesPyramid::E_E ( int p,
Vector s,
Vector sds,
DenseMatrix & u ) const

This is a vector-valued function associated with an edge of a pyramid

The vector s contains two coordinate values and ds is related to the gradient of these coordinates with respect to the reference coordinates i.e. sds = s0 grad s1 - s1 grad s0

Definition at line 960 of file fe_pyramid.cpp.

◆ E_Q() [1/2]

void mfem::FuentesPyramid::E_Q ( int p,
Vector s,
const DenseMatrix & grad_s,
Vector t,
const DenseMatrix & grad_t,
DenseTensor & u,
DenseTensor & curl_u ) const

Definition at line 1066 of file fe_pyramid.cpp.

◆ E_Q() [2/2]

void mfem::FuentesPyramid::E_Q ( int p,
Vector s,
Vector ds,
Vector t,
DenseTensor & u ) const

Definition at line 1027 of file fe_pyramid.cpp.

◆ E_T() [1/2]

void mfem::FuentesPyramid::E_T ( int p,
Vector s,
const DenseMatrix & grad_s,
DenseTensor & u,
DenseTensor & curl_u ) const

Definition at line 1171 of file fe_pyramid.cpp.

◆ E_T() [2/2]

void mfem::FuentesPyramid::E_T ( int p,
Vector s,
Vector sds,
DenseTensor & u ) const

Definition at line 1137 of file fe_pyramid.cpp.

◆ grad_lam1()

Vector mfem::FuentesPyramid::grad_lam1 ( real_t x,
real_t y,
real_t z )
static

Gradients of the "Affine" Coordinates.

Definition at line 21 of file fe_pyramid.cpp.

◆ grad_lam15()

DenseMatrix mfem::FuentesPyramid::grad_lam15 ( real_t x,
real_t y,
real_t z )
static

Gradients of the above two component vectors.

Definition at line 54 of file fe_pyramid.cpp.

◆ grad_lam2()

Vector mfem::FuentesPyramid::grad_lam2 ( real_t x,
real_t y,
real_t z )
static

Definition at line 28 of file fe_pyramid.cpp.

◆ grad_lam25()

DenseMatrix mfem::FuentesPyramid::grad_lam25 ( real_t x,
real_t y,
real_t z )
static

Definition at line 62 of file fe_pyramid.cpp.

◆ grad_lam3()

Vector mfem::FuentesPyramid::grad_lam3 ( real_t x,
real_t y,
real_t z )
static

Definition at line 35 of file fe_pyramid.cpp.

◆ grad_lam35()

DenseMatrix mfem::FuentesPyramid::grad_lam35 ( real_t x,
real_t y,
real_t z )
static

Definition at line 70 of file fe_pyramid.cpp.

◆ grad_lam4()

Vector mfem::FuentesPyramid::grad_lam4 ( real_t x,
real_t y,
real_t z )
static

Definition at line 42 of file fe_pyramid.cpp.

◆ grad_lam45()

DenseMatrix mfem::FuentesPyramid::grad_lam45 ( real_t x,
real_t y,
real_t z )
static

Definition at line 78 of file fe_pyramid.cpp.

◆ grad_lam5()

Vector mfem::FuentesPyramid::grad_lam5 ( real_t x,
real_t y,
real_t z )
static

Definition at line 49 of file fe_pyramid.cpp.

◆ grad_mu0() [1/2]

static Vector mfem::FuentesPyramid::grad_mu0 ( real_t z)
inlinestatic

Definition at line 172 of file fe_pyramid.hpp.

◆ grad_mu0() [2/2]

Vector mfem::FuentesPyramid::grad_mu0 ( real_t z,
const Vector xy,
unsigned int ab )
static

Definition at line 186 of file fe_pyramid.cpp.

◆ grad_mu01() [1/2]

DenseMatrix mfem::FuentesPyramid::grad_mu01 ( real_t z)
static

Definition at line 178 of file fe_pyramid.cpp.

◆ grad_mu01() [2/2]

DenseMatrix mfem::FuentesPyramid::grad_mu01 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 200 of file fe_pyramid.cpp.

◆ grad_mu1() [1/2]

static Vector mfem::FuentesPyramid::grad_mu1 ( real_t z)
inlinestatic

Definition at line 174 of file fe_pyramid.hpp.

◆ grad_mu1() [2/2]

Vector mfem::FuentesPyramid::grad_mu1 ( real_t z,
const Vector xy,
unsigned int ab )
static

Definition at line 193 of file fe_pyramid.cpp.

◆ grad_nu0()

Vector mfem::FuentesPyramid::grad_nu0 ( real_t z,
const Vector xy,
unsigned int ab )
static

Definition at line 216 of file fe_pyramid.cpp.

◆ grad_nu01()

DenseMatrix mfem::FuentesPyramid::grad_nu01 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 233 of file fe_pyramid.cpp.

◆ grad_nu012()

DenseMatrix mfem::FuentesPyramid::grad_nu012 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 241 of file fe_pyramid.cpp.

◆ grad_nu1()

Vector mfem::FuentesPyramid::grad_nu1 ( real_t z,
const Vector xy,
unsigned int ab )
static

Definition at line 222 of file fe_pyramid.cpp.

◆ grad_nu120()

DenseMatrix mfem::FuentesPyramid::grad_nu120 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 250 of file fe_pyramid.cpp.

◆ grad_nu2()

Vector mfem::FuentesPyramid::grad_nu2 ( real_t z,
const Vector xy,
unsigned int ab )
static

Definition at line 228 of file fe_pyramid.cpp.

◆ lam1()

static real_t mfem::FuentesPyramid::lam1 ( real_t x,
real_t y,
real_t z )
inlinestatic

Pyramid "Affine" Coordinates.

Definition at line 90 of file fe_pyramid.hpp.

◆ lam125()

static Vector mfem::FuentesPyramid::lam125 ( real_t x,
real_t y,
real_t z )
inlinestatic

Three component vectors associated with triangular faces.

Definition at line 131 of file fe_pyramid.hpp.

◆ lam125_grad_lam125()

Vector mfem::FuentesPyramid::lam125_grad_lam125 ( real_t x,
real_t y,
real_t z )
static

Vector functions related to the normals to the triangular faces

Computes λiλj×λ5+λjλ5×λi+λ5λi×λj

Definition at line 118 of file fe_pyramid.cpp.

◆ lam145()

static Vector mfem::FuentesPyramid::lam145 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 141 of file fe_pyramid.hpp.

◆ lam145_grad_lam145()

Vector mfem::FuentesPyramid::lam145_grad_lam145 ( real_t x,
real_t y,
real_t z )
static

Definition at line 153 of file fe_pyramid.cpp.

◆ lam15()

static Vector mfem::FuentesPyramid::lam15 ( real_t x,
real_t y,
real_t z )
inlinestatic

Two component vectors associated with edges touching the apex.

Definition at line 109 of file fe_pyramid.hpp.

◆ lam15_grad_lam15()

Vector mfem::FuentesPyramid::lam15_grad_lam15 ( real_t x,
real_t y,
real_t z )
static

Computes λiλ5λ5λi.

Definition at line 86 of file fe_pyramid.cpp.

◆ lam2()

static real_t mfem::FuentesPyramid::lam2 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 92 of file fe_pyramid.hpp.

◆ lam235()

static Vector mfem::FuentesPyramid::lam235 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 133 of file fe_pyramid.hpp.

◆ lam235_grad_lam235()

Vector mfem::FuentesPyramid::lam235_grad_lam235 ( real_t x,
real_t y,
real_t z )
static

Definition at line 125 of file fe_pyramid.cpp.

◆ lam25()

static Vector mfem::FuentesPyramid::lam25 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 111 of file fe_pyramid.hpp.

◆ lam25_grad_lam25()

Vector mfem::FuentesPyramid::lam25_grad_lam25 ( real_t x,
real_t y,
real_t z )
static

Definition at line 94 of file fe_pyramid.cpp.

◆ lam3()

static real_t mfem::FuentesPyramid::lam3 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 94 of file fe_pyramid.hpp.

◆ lam345()

static Vector mfem::FuentesPyramid::lam345 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 135 of file fe_pyramid.hpp.

◆ lam345_grad_lam345()

Vector mfem::FuentesPyramid::lam345_grad_lam345 ( real_t x,
real_t y,
real_t z )
static

Definition at line 132 of file fe_pyramid.cpp.

◆ lam35()

static Vector mfem::FuentesPyramid::lam35 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 113 of file fe_pyramid.hpp.

◆ lam35_grad_lam35()

Vector mfem::FuentesPyramid::lam35_grad_lam35 ( real_t x,
real_t y,
real_t z )
static

Definition at line 102 of file fe_pyramid.cpp.

◆ lam4()

static real_t mfem::FuentesPyramid::lam4 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 96 of file fe_pyramid.hpp.

◆ lam415()

static Vector mfem::FuentesPyramid::lam415 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 139 of file fe_pyramid.hpp.

◆ lam415_grad_lam415()

Vector mfem::FuentesPyramid::lam415_grad_lam415 ( real_t x,
real_t y,
real_t z )
static

Definition at line 146 of file fe_pyramid.cpp.

◆ lam435()

static Vector mfem::FuentesPyramid::lam435 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 137 of file fe_pyramid.hpp.

◆ lam435_grad_lam435()

Vector mfem::FuentesPyramid::lam435_grad_lam435 ( real_t x,
real_t y,
real_t z )
static

Definition at line 139 of file fe_pyramid.cpp.

◆ lam45()

static Vector mfem::FuentesPyramid::lam45 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 115 of file fe_pyramid.hpp.

◆ lam45_grad_lam45()

Vector mfem::FuentesPyramid::lam45_grad_lam45 ( real_t x,
real_t y,
real_t z )
static

Definition at line 110 of file fe_pyramid.cpp.

◆ lam5()

static real_t mfem::FuentesPyramid::lam5 ( real_t x,
real_t y,
real_t z )
inlinestatic

Definition at line 98 of file fe_pyramid.hpp.

◆ mu0() [1/2]

static real_t mfem::FuentesPyramid::mu0 ( real_t z)
inlinestatic

Definition at line 167 of file fe_pyramid.hpp.

◆ mu0() [2/2]

static real_t mfem::FuentesPyramid::mu0 ( real_t z,
const Vector & xy,
unsigned int ab )
inlinestatic

Definition at line 182 of file fe_pyramid.hpp.

◆ mu01() [1/2]

static Vector mfem::FuentesPyramid::mu01 ( real_t z)
inlinestatic

Definition at line 177 of file fe_pyramid.hpp.

◆ mu01() [2/2]

static Vector mfem::FuentesPyramid::mu01 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 190 of file fe_pyramid.hpp.

◆ mu01_grad_mu01()

Vector mfem::FuentesPyramid::mu01_grad_mu01 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 208 of file fe_pyramid.cpp.

◆ mu1() [1/2]

static real_t mfem::FuentesPyramid::mu1 ( real_t z)
inlinestatic

Definition at line 169 of file fe_pyramid.hpp.

◆ mu1() [2/2]

static real_t mfem::FuentesPyramid::mu1 ( real_t z,
const Vector & xy,
unsigned int ab )
inlinestatic

Definition at line 184 of file fe_pyramid.hpp.

◆ nu0()

static real_t mfem::FuentesPyramid::nu0 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 196 of file fe_pyramid.hpp.

◆ nu01()

static Vector mfem::FuentesPyramid::nu01 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 205 of file fe_pyramid.hpp.

◆ nu012()

static Vector mfem::FuentesPyramid::nu012 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 209 of file fe_pyramid.hpp.

◆ nu012_grad_nu012()

Vector mfem::FuentesPyramid::nu012_grad_nu012 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 275 of file fe_pyramid.cpp.

◆ nu01_grad_nu01()

Vector mfem::FuentesPyramid::nu01_grad_nu01 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 259 of file fe_pyramid.cpp.

◆ nu1()

static real_t mfem::FuentesPyramid::nu1 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 198 of file fe_pyramid.hpp.

◆ nu12()

static Vector mfem::FuentesPyramid::nu12 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 207 of file fe_pyramid.hpp.

◆ nu120()

static Vector mfem::FuentesPyramid::nu120 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 211 of file fe_pyramid.hpp.

◆ nu12_grad_nu12()

Vector mfem::FuentesPyramid::nu12_grad_nu12 ( real_t z,
Vector xy,
unsigned int ab )
static

Definition at line 267 of file fe_pyramid.cpp.

◆ nu2()

static real_t mfem::FuentesPyramid::nu2 ( real_t z,
Vector xy,
unsigned int ab )
inlinestatic

Definition at line 199 of file fe_pyramid.hpp.

◆ phi_E() [1/5]

void mfem::FuentesPyramid::phi_E ( int p,
real_t s0,
real_t s1,
real_t * u )
static

u must be at least p + 1 in length

Definition at line 730 of file fe_pyramid.cpp.

◆ phi_E() [2/5]

void mfem::FuentesPyramid::phi_E ( int p,
real_t s0,
real_t s1,
real_t * u,
real_t * duds0,
real_t * duds1 )
static

Definition at line 737 of file fe_pyramid.cpp.

◆ phi_E() [3/5]

void mfem::FuentesPyramid::phi_E ( int p,
Vector s,
const DenseMatrix & grad_s,
Vector & u,
DenseMatrix & grad_u ) const

grad_s must be 2x3

Definition at line 762 of file fe_pyramid.cpp.

◆ phi_E() [4/5]

void mfem::FuentesPyramid::phi_E ( int p,
Vector s,
Vector & u )
static

Definition at line 744 of file fe_pyramid.cpp.

◆ phi_E() [5/5]

void mfem::FuentesPyramid::phi_E ( int p,
Vector s,
Vector & u,
DenseMatrix & duds )
static

Definition at line 751 of file fe_pyramid.cpp.

◆ phi_Q() [1/2]

void mfem::FuentesPyramid::phi_Q ( int p,
Vector s,
const DenseMatrix & grad_s,
Vector t,
const DenseMatrix & grad_t,
DenseMatrix & u,
DenseTensor & grad_u ) const

Definition at line 814 of file fe_pyramid.cpp.

◆ phi_Q() [2/2]

void mfem::FuentesPyramid::phi_Q ( int p,
Vector s,
Vector t,
DenseMatrix & u ) const

u must be at least (p+1)x(p+1) in size

Definition at line 786 of file fe_pyramid.cpp.

◆ phi_T() [1/2]

void mfem::FuentesPyramid::phi_T ( int p,
Vector nu,
const DenseMatrix & grad_nu,
DenseMatrix & u,
DenseTensor & grad_u ) const

Definition at line 900 of file fe_pyramid.cpp.

◆ phi_T() [2/2]

void mfem::FuentesPyramid::phi_T ( int p,
Vector nu,
DenseMatrix & u ) const

Definition at line 868 of file fe_pyramid.cpp.

◆ V_L()

void mfem::FuentesPyramid::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

This implements Vij from the Fuentes paper

u must be at least (p+1)x(p+1)x3

Definition at line 1469 of file fe_pyramid.cpp.

◆ V_Q()

void mfem::FuentesPyramid::V_Q ( int p,
Vector s,
Vector ds,
Vector t,
Vector dt,
DenseTensor & u ) const

This is a vector-valued function associated with the quadrilateral face of a pyramid

The vectors s and t contain pairs of coordinate values and ds and dt are related to derivatives of these coordinates:

ds = s0 grad s1 - s1 grad s0

dt = t0 grad t1 - t1 grad t0

Definition at line 1238 of file fe_pyramid.cpp.

◆ V_R()

void mfem::FuentesPyramid::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

This implements Vi from the Fuentes paper

u must be at least (p+1)x3

Definition at line 1538 of file fe_pyramid.cpp.

◆ V_T() [1/2]

void mfem::FuentesPyramid::V_T ( int p,
Vector s,
Vector sdsxds,
DenseTensor & u ) const

This is a vector-valued function associated with the triangular faces of a pyramid

The vector s contains three coordinate values and sdsxds is related to derivatives of these coordinates with respect to the reference coordinates:

sdsxds = s0 grad s1 x grad s2 + s1 grad s2 x grad s0 + s2 grad s0 x grad s1

Definition at line 1272 of file fe_pyramid.cpp.

◆ V_T() [2/2]

void mfem::FuentesPyramid::V_T ( int p,
Vector s,
Vector sdsxds,
real_t dsdsxds,
DenseTensor & u,
DenseMatrix & du ) const

This computes V_T as above and its divergence

The vector s contains three coordinate values and sdsxds is related to derivatives of these coordinates with respect to the reference coordinates:

sdsxds = s0 grad s1 x grad s2 + s1 grad s2 x grad s0 + s2 grad s0 x grad s1

The scalar dsdsxds is the divergence of sdsxds:

dsdsxds = grad s0 dot (grad s1 x grad s2)

Definition at line 1306 of file fe_pyramid.cpp.

◆ VT_T() [1/2]

void mfem::FuentesPyramid::VT_T ( int p,
Vector s,
Vector sds,
Vector sdsxds,
real_t mu,
Vector grad_mu,
DenseTensor & u ) const

Definition at line 1345 of file fe_pyramid.cpp.

◆ VT_T() [2/2]

void mfem::FuentesPyramid::VT_T ( int p,
Vector s,
Vector sds,
Vector sdsxds,
Vector grad_s2,
real_t mu,
Vector grad_mu,
DenseTensor & u,
DenseMatrix & du ) const

Definition at line 1398 of file fe_pyramid.cpp.

Member Data Documentation

◆ apex_tol

real_t mfem::FuentesPyramid::apex_tol = 1e-8
staticconstexprprotected

Definition at line 82 of file fe_pyramid.hpp.

◆ one

real_t mfem::FuentesPyramid::one = 1.0
staticconstexprprotected

Definition at line 80 of file fe_pyramid.hpp.

◆ zero

real_t mfem::FuentesPyramid::zero = 0.0
staticconstexprprotected

Definition at line 81 of file fe_pyramid.hpp.


The documentation for this class was generated from the following files: