Loading...
Searching...
No Matches
terra::fe::wedge Namespace Reference

Features for wedge elements. More...

Namespaces

namespace  linearforms
 
namespace  operators
 
namespace  quadrature
 Quadrature rules for the reference wedge.
 
namespace  shell
 

Enumerations

enum class  DoFOrdering { NODEWISE , DIMENSIONWISE }
 

Functions

template<std::floating_point T>
constexpr T shape_rad (const int node_idx, const T zeta)
 Radial shape function.
 
template<std::floating_point T>
constexpr T shape_rad (const int node_idx, const dense::Vec< T, 3 > &xi_eta_zeta)
 Radial shape function.
 
template<std::floating_point T>
constexpr T shape_lat (const int node_idx, const T xi, const T eta)
 Lateral shape function.
 
template<std::floating_point T>
constexpr T shape_lat (const int node_idx, const dense::Vec< T, 3 > &xi_eta_zeta)
 Lateral shape function.
 
template<std::floating_point T>
constexpr T shape (const int node_idx, const T xi, const T eta, const T zeta)
 (Tensor-product) Shape function.
 
template<std::floating_point T>
constexpr T shape (const int node_idx, const dense::Vec< T, 3 > &xi_eta_zeta)
 (Tensor-product) Shape function.
 
template<std::floating_point T>
constexpr T grad_shape_rad (const int node_idx)
 Gradient of the radial part of the shape function, in the radial direction \( \frac{\partial}{\partial \zeta} N^\mathrm{rad}_j \).
 
template<std::floating_point T>
constexpr T grad_shape_lat_xi (const int node_idx)
 Gradient of the lateral part of the shape function, in xi direction \( \frac{\partial}{\partial \xi} N^\mathrm{lat}_j \).
 
template<std::floating_point T>
constexpr T grad_shape_lat_eta (const int node_idx)
 Gradient of the lateral part of the shape function, in eta direction \( \frac{\partial}{\partial \eta} N^\mathrm{lat}_j \).
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > grad_shape (const int node_idx, const T xi, const T eta, const T zeta)
 Gradient of the full shape function:
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > grad_shape (const int node_idx, const dense::Vec< T, 3 > &xi_eta_zeta)
 Gradient of the full shape function:
 
template<std::floating_point T>
constexpr T shape_rad_coarse (const int coarse_node_idx, const int fine_radial_wedge_idx, const T zeta_fine)
 Returns the coarse grid radial shape function evaluated at a point of the reference fine grid wedge.
 
template<std::floating_point T>
constexpr T shape_lat_coarse (const int coarse_node_idx, const int fine_lateral_wedge_idx, const T xi_fine, const T eta_fine)
 Returns the coarse grid lateral shape function evaluated at a point of the reference fine grid wedge.
 
template<std::floating_point T>
constexpr T shape_coarse (const int coarse_node_idx, const int fine_radial_wedge_idx, const int fine_lateral_wedge_idx, const T xi_fine, const T eta_fine, const T zeta_fine)
 
template<std::floating_point T>
constexpr T shape_coarse (const int coarse_node_idx, const int fine_radial_wedge_idx, const int fine_lateral_wedge_idx, const dense::Vec< T, 3 > &xi_eta_zeta_fine)
 
template<std::floating_point T>
constexpr T grad_shape_rad_coarse (const int coarse_node_idx, const int fine_radial_wedge_idx)
 
template<std::floating_point T>
constexpr T grad_shape_lat_coarse_xi (const int coarse_node_idx, const int fine_lateral_wedge_idx)
 
template<std::floating_point T>
constexpr T grad_shape_lat_coarse_eta (const int coarse_node_idx, const int fine_lateral_wedge_idx)
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > grad_shape_coarse (const int node_idx, const int fine_radial_wedge_idx, const int fine_lateral_wedge_idx, const T xi, const T eta, const T zeta)
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > grad_shape_coarse (const int node_idx, const int fine_radial_wedge_idx, const int fine_lateral_wedge_idx, const dense::Vec< T, 3 > &xi_eta_zeta_fine)
 
template<std::floating_point T>
constexpr T forward_map_rad (const T r_1, const T r_2, const T zeta)
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > forward_map_lat (const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &p2_phy, const dense::Vec< T, 3 > &p3_phy, const T xi, const T eta)
 
template<std::floating_point T>
constexpr T grad_forward_map_rad (const T r_1, const T r_2)
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > grad_forward_map_lat_xi (const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &p2_phy, const dense::Vec< T, 3 > &)
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > grad_forward_map_lat_eta (const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &, const dense::Vec< T, 3 > &p3_phy)
 
template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > jac_lat (const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &p2_phy, const dense::Vec< T, 3 > &p3_phy, const T xi, const T eta)
 
template<std::floating_point T>
constexpr dense::Vec< T, 3 > jac_rad (const T r_1, const T r_2, const T zeta)
 
template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > jac (const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &p2_phy, const dense::Vec< T, 3 > &p3_phy, const T r_1, const T r_2, const T xi, const T eta, const T zeta)
 
template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > jac (const dense::Vec< T, 3 > p_phy[3], const T r_1, const T r_2, const dense::Vec< T, 3 > &xi_eta_zeta_fine)
 
template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > symmetric_grad (const dense::Mat< T, 3, 3 > &J_inv_transposed, const dense::Vec< T, 3 > &quad_point, const int dof, const int dim)
 Returns the symmetric gradient of the shape function of a dof at a quadrature point.
 
template<std::floating_point T>
void wedge_surface_physical_coords (dense::Vec< T, 3 >(&wedge_surf_phy_coords)[num_wedges_per_hex_cell][num_nodes_per_wedge_surface], const grid::Grid3DDataVec< T, 3 > &lateral_grid, const int local_subdomain_id, const int x_cell, const int y_cell)
 Extracts the (unit sphere) surface vertex coords of the two wedges of a hex cell.
 
template<std::floating_point T>
void wedge_0_surface_physical_coords (dense::Vec< T, 3 > *wedge_surf_phy_coords, const grid::Grid3DDataVec< T, 3 > &lateral_grid, const int local_subdomain_id, const int x_cell, const int y_cell)
 
template<std::floating_point T>
void wedge_1_surface_physical_coords (dense::Vec< T, 3 > *wedge_surf_phy_coords, const grid::Grid3DDataVec< T, 3 > &lateral_grid, const int local_subdomain_id, const int x_cell, const int y_cell)
 
template<std::floating_point T, int NumQuadPoints>
constexpr void jacobian_lat_inverse_transposed_and_determinant (dense::Mat< T, 3, 3 >(&jac_lat_inv_t)[num_wedges_per_hex_cell][NumQuadPoints], T(&det_jac_lat)[num_wedges_per_hex_cell][NumQuadPoints], const dense::Vec< T, 3 > wedge_surf_phy_coords[num_wedges_per_hex_cell][num_nodes_per_wedge_surface], const dense::Vec< T, 3 > quad_points[NumQuadPoints])
 Computes the transposed inverse of the Jacobian of the lateral forward map from the reference triangle to the triangle on the unit sphere and the absolute determinant of that Jacobian at the passed quadrature points.
 
template<std::floating_point T, int NumQuadPoints>
constexpr void jacobian_lat_determinant (T(&det_jac_lat)[num_wedges_per_hex_cell][NumQuadPoints], const dense::Vec< T, 3 > wedge_surf_phy_coords[num_wedges_per_hex_cell][num_nodes_per_wedge_surface], const dense::Vec< T, 3 > quad_points[NumQuadPoints])
 Like jacobian_lat_inverse_transposed_and_determinant() but only computes the determinant (cheaper if the transposed inverse of the Jacobian is not required).
 
template<std::floating_point T, int NumQuadPoints>
constexpr void lateral_parts_of_grad_phi (dense::Vec< T, 3 >(&g_rad)[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints], dense::Vec< T, 3 >(&g_lat)[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints], const dense::Mat< T, 3, 3 > jac_lat_inv_t[num_wedges_per_hex_cell][NumQuadPoints], const dense::Vec< T, 3 > quad_points[NumQuadPoints])
 Computes the radially independent parts of the physical shape function gradients.
 
template<std::floating_point T, int NumQuadPoints>
constexpr dense::Vec< T, 3 > grad_shape_full (const dense::Vec< T, 3 > g_rad[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints], const dense::Vec< T, 3 > g_lat[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints], const T r_inv, const T grad_r_inv, const int wedge_idx, const int node_idx, const int quad_point_idx)
 Computes and returns J^-T grad(N_j).
 
template<std::floating_point T, int NumQuadPoints>
constexpr T det_full (const T det_jac_lat[num_wedges_per_hex_cell][NumQuadPoints], const T r, const T grad_r, const int wedge_idx, const int quad_point_idx)
 Computes |det(J)|.
 
template<std::floating_point T>
void extract_local_wedge_scalar_coefficients (dense::Vec< T, 6 >(&local_coefficients)[2], const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const grid::Grid4DDataScalar< T > &global_coefficients)
 Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient vector.
 
template<std::floating_point T, int VecDim>
void extract_local_wedge_vector_coefficients (dense::Vec< T, 6 >(&local_coefficients)[2], const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int d, const grid::Grid4DDataVec< T, VecDim > &global_coefficients)
 Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient vector.
 
template<std::floating_point T>
void atomically_add_local_wedge_scalar_coefficients (const grid::Grid4DDataScalar< T > &global_coefficients, const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const dense::Vec< T, 6 >(&local_coefficients)[2])
 Performs an atomic add of the two local wedge coefficient vectors of a hex cell into the global coefficient vector.
 
template<std::floating_point T, int VecDim>
void atomically_add_local_wedge_vector_coefficients (const grid::Grid4DDataVec< T, VecDim > &global_coefficients, const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int d, const dense::Vec< T, 6 > local_coefficients[2])
 Performs an atomic add of the two local wedge coefficient vectors of a hex cell into the global coefficient vector.
 
constexpr int fine_lateral_wedge_idx (const int x_cell_fine, const int y_cell_fine, const int wedge_idx_fine)
 Returns the lateral wedge index with respect to a coarse grid wedge from the fine wedge indices.
 
template<typename ScalarT >
constexpr void reorder_local_dofs (const DoFOrdering doo_from, const DoFOrdering doo_to, dense::Vec< ScalarT, 18 > &dofs)
 

Variables

constexpr int num_wedges_per_hex_cell = 2
 
constexpr int num_nodes_per_wedge_surface = 3
 
constexpr int num_nodes_per_wedge = 6
 

Detailed Description

Features for wedge elements.

Geometry

\( \xi, \eta \in [0, 1] \) (lateral reference coords)

\( \zeta \in [-1, 1] \) (radial reference coords)

\( 0 \leq \xi + \eta \leq 1 \)

Node enumeration

r_node_idx = r_cell_idx + 1 (outer):
5
|\
| \
3--4
r_node_idx = r_cell_idx (inner):
2
|\
| \
0--1

Shape functions

Lateral:

\[ \begin{align} N^\mathrm{lat}_0 = N^\mathrm{lat}_3 &= 1 - \xi - \eta \\ N^\mathrm{lat}_1 = N^\mathrm{lat}_4 &= \xi \\ N^\mathrm{lat}_2 = N^\mathrm{lat}_5 &= \eta \end{align} \]

Radial:

\[ \begin{align} N^\mathrm{rad}_0 = N^\mathrm{rad}_1 = N^\mathrm{rad}_2 &= \frac{1}{2} ( 1 - \zeta ) \\ N^\mathrm{rad}_3 = N^\mathrm{rad}_4 = N^\mathrm{rad}_5 &= \frac{1}{2} ( 1 + \zeta ) \\ \end{align} \]

Full:

\[ N_i = N^\mathrm{lat}_i * N^\mathrm{rad}_i \]

Physical coordinates

r_1, r_2 radii of bottom and top (r_1 < r_2)
p1_phy, p2_phy, p3_phy coords of triangle on unit sphere

Enumeration Type Documentation

◆ DoFOrdering

enum class terra::fe::wedge::DoFOrdering
strong
Enumerator
NODEWISE 
DIMENSIONWISE 

Function Documentation

◆ atomically_add_local_wedge_scalar_coefficients()

template<std::floating_point T>
void terra::fe::wedge::atomically_add_local_wedge_scalar_coefficients ( const grid::Grid4DDataScalar< T > &  global_coefficients,
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell,
const int  r_cell,
const dense::Vec< T, 6 >(&)  local_coefficients[2] 
)

Performs an atomic add of the two local wedge coefficient vectors of a hex cell into the global coefficient vector.

r = r_cell + 1 (outer)
6--7
|\ |
| \|
4--5
r = r_cell (inner)
2--3
|\ |
| \|
0--1
v0 = (0, 1, 2, 4, 5, 6)
v1 = (3, 2, 1, 7, 6, 5)
Parameters
global_coefficients[inout] the global coefficient vector
local_subdomain_id[in] shell subdomain id on this process
x_cell[in] hex cell x-coordinate
y_cell[in] hex cell y-coordinate
r_cell[in] hex cell r-coordinate
local_coefficients[in] the local coefficient vector

◆ atomically_add_local_wedge_vector_coefficients()

template<std::floating_point T, int VecDim>
void terra::fe::wedge::atomically_add_local_wedge_vector_coefficients ( const grid::Grid4DDataVec< T, VecDim > &  global_coefficients,
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell,
const int  r_cell,
const int  d,
const dense::Vec< T, 6 >  local_coefficients[2] 
)

Performs an atomic add of the two local wedge coefficient vectors of a hex cell into the global coefficient vector.

r = r_cell + 1 (outer)
6--7
|\ |
| \|
4--5
r = r_cell (inner)
2--3
|\ |
| \|
0--1
v0 = (0, 1, 2, 4, 5, 6)
v1 = (3, 2, 1, 7, 6, 5)
Parameters
global_coefficients[inout] the global coefficient vector
local_subdomain_id[in] shell subdomain id on this process
x_cell[in] hex cell x-coordinate
y_cell[in] hex cell y-coordinate
r_cell[in] hex cell r-coordinate
d[in] vector-element of the vector-valued global view
local_coefficients[in] the local coefficient vector

◆ det_full()

template<std::floating_point T, int NumQuadPoints>
constexpr T terra::fe::wedge::det_full ( const T  det_jac_lat[num_wedges_per_hex_cell][NumQuadPoints],
const T  r,
const T  grad_r,
const int  wedge_idx,
const int  quad_point_idx 
)
constexpr

Computes |det(J)|.

Parameters
det_jac_lat[in] see jacobian_lat_determinant()
r[in] the physical radius of the quadrature point - compute with forward_map_rad()
grad_r[in] the radial component of the gradient of the forward map in radial direction - compute with grad_forward_map_rad())
wedge_idx[in] wedge index of the hex cell (0 or 1)
quad_point_idx[in] index of the quadrature point

◆ extract_local_wedge_scalar_coefficients()

template<std::floating_point T>
void terra::fe::wedge::extract_local_wedge_scalar_coefficients ( dense::Vec< T, 6 >(&)  local_coefficients[2],
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell,
const int  r_cell,
const grid::Grid4DDataScalar< T > &  global_coefficients 
)

Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient vector.

r = r_cell + 1 (outer)
6--7
|\ |
| \|
4--5
r = r_cell (inner)
2--3
|\ |
| \|
0--1
v0 = (0, 1, 2, 4, 5, 6)
v1 = (3, 2, 1, 7, 6, 5)
Parameters
local_coefficients[out] the local coefficient vector
local_subdomain_id[in] shell subdomain id on this process
x_cell[in] hex cell x-coordinate
y_cell[in] hex cell y-coordinate
r_cell[in] hex cell r-coordinate
global_coefficients[in] the global coefficient vector

◆ extract_local_wedge_vector_coefficients()

template<std::floating_point T, int VecDim>
void terra::fe::wedge::extract_local_wedge_vector_coefficients ( dense::Vec< T, 6 >(&)  local_coefficients[2],
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell,
const int  r_cell,
const int  d,
const grid::Grid4DDataVec< T, VecDim > &  global_coefficients 
)

Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient vector.

r = r_cell + 1 (outer)
6--7
|\ |
| \|
4--5
r = r_cell (inner)
2--3
|\ |
| \|
0--1
v0 = (0, 1, 2, 4, 5, 6)
v1 = (3, 2, 1, 7, 6, 5)
Parameters
local_coefficients[out] the local coefficient vector
local_subdomain_id[in] shell subdomain id on this process
x_cell[in] hex cell x-coordinate
y_cell[in] hex cell y-coordinate
r_cell[in] hex cell r-coordinate
d[in] vector-element of the vector-valued global view
global_coefficients[in] the global coefficient vector

◆ fine_lateral_wedge_idx()

constexpr int terra::fe::wedge::fine_lateral_wedge_idx ( const int  x_cell_fine,
const int  y_cell_fine,
const int  wedge_idx_fine 
)
constexpr

Returns the lateral wedge index with respect to a coarse grid wedge from the fine wedge indices.

Each coarse grid wedge is laterally divided into four fine wedges (radially into two). The lateral four fine wedges are enumerated from 0 to 3.

This function returns that lateral index given a fine grid index of a wedge.

Here are two coarse wedges (view from the "top"), each with four fine grid wedges (enumerated from 0 to 3).

Coarse wedge idx = 0
+
|\
| \
| \
| 2 \
+----+
|\ 3 |\
| \ | \
| \ | \
| 0 \| 1 \
+----+----+
Coarse wedge idx = 1
+----+----+
\ 1 |\ 0 |
\ | \ |
\ | \ |
\| 3 \|
+----+
\ 2 |
\ |
\ |
\|
+

This function now computes the indices from a grid of fine wedges that looks like this:

x_cell_fine % 1:
+----+----+
|\ 0 |\ 1 |
| \ | \ |
| \ | \ |
| 0 \| 1 \|
+----+----+
|\ 0 |\ 1 |
| \ | \ |
| \ | \ |
| 0 \| 1 \|
+----+----+
y_cell_fine % 1:
+----+----+
|\ 1 |\ 1 |
| \ | \ |
| \ | \ |
| 1 \| 1 \|
+----+----+
|\ 0 |\ 0 |
| \ | \ |
| \ | \ |
| 0 \| 0 \|
+----+----+
wedge_idx_fine:
+----+----+
|\ 1 |\ 1 |
| \ | \ |
| \ | \ |
| 0 \| 0 \|
+----+----+
|\ 1 |\ 1 |
| \ | \ |
| \ | \ |
| 0 \| 0 \|
+----+----+
Resulting/returned values:
+----+----+
|\ 1 |\ 0 |
| \ | \ |
| \ | \ |
| 2 \| 3 \|
+----+----+
|\ 3 |\ 2 |
| \ | \ |
| \ | \ |
| 0 \| 1 \|
+----+----+

◆ forward_map_lat()

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::forward_map_lat ( const dense::Vec< T, 3 > &  p1_phy,
const dense::Vec< T, 3 > &  p2_phy,
const dense::Vec< T, 3 > &  p3_phy,
const T  xi,
const T  eta 
)
constexpr

◆ forward_map_rad()

template<std::floating_point T>
constexpr T terra::fe::wedge::forward_map_rad ( const T  r_1,
const T  r_2,
const T  zeta 
)
constexpr

◆ grad_forward_map_lat_eta()

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_forward_map_lat_eta ( const dense::Vec< T, 3 > &  p1_phy,
const dense::Vec< T, 3 > &  ,
const dense::Vec< T, 3 > &  p3_phy 
)
constexpr

◆ grad_forward_map_lat_xi()

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_forward_map_lat_xi ( const dense::Vec< T, 3 > &  p1_phy,
const dense::Vec< T, 3 > &  p2_phy,
const dense::Vec< T, 3 > &   
)
constexpr

◆ grad_forward_map_rad()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_forward_map_rad ( const T  r_1,
const T  r_2 
)
constexpr

◆ grad_shape() [1/2]

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_shape ( const int  node_idx,
const dense::Vec< T, 3 > &  xi_eta_zeta 
)
constexpr

Gradient of the full shape function:

\[ \nabla N_j = \begin{bmatrix} \frac{\partial}{\partial \xi} N_j \\ \frac{\partial}{\partial \eta} N_j \\ \frac{\partial}{\partial \zeta} N_j \end{bmatrix} = \begin{bmatrix} N^\mathrm{rad}_j \frac{\partial}{\partial \xi} N^\mathrm{lat}_j \\ N^\mathrm{rad}_j \frac{\partial}{\partial \eta} N^\mathrm{lat}_j \\ N^\mathrm{lat}_j \frac{\partial}{\partial \zeta} N^\mathrm{rad}_j \end{bmatrix} \]

◆ grad_shape() [2/2]

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_shape ( const int  node_idx,
const T  xi,
const T  eta,
const T  zeta 
)
constexpr

Gradient of the full shape function:

\[ \nabla N_j = \begin{bmatrix} \frac{\partial}{\partial \xi} N_j \\ \frac{\partial}{\partial \eta} N_j \\ \frac{\partial}{\partial \zeta} N_j \end{bmatrix} = \begin{bmatrix} N^\mathrm{rad}_j \frac{\partial}{\partial \xi} N^\mathrm{lat}_j \\ N^\mathrm{rad}_j \frac{\partial}{\partial \eta} N^\mathrm{lat}_j \\ N^\mathrm{lat}_j \frac{\partial}{\partial \zeta} N^\mathrm{rad}_j \end{bmatrix} \]

◆ grad_shape_coarse() [1/2]

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_shape_coarse ( const int  node_idx,
const int  fine_radial_wedge_idx,
const int  fine_lateral_wedge_idx,
const dense::Vec< T, 3 > &  xi_eta_zeta_fine 
)
constexpr

◆ grad_shape_coarse() [2/2]

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_shape_coarse ( const int  node_idx,
const int  fine_radial_wedge_idx,
const int  fine_lateral_wedge_idx,
const T  xi,
const T  eta,
const T  zeta 
)
constexpr

◆ grad_shape_full()

template<std::floating_point T, int NumQuadPoints>
constexpr dense::Vec< T, 3 > terra::fe::wedge::grad_shape_full ( const dense::Vec< T, 3 >  g_rad[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints],
const dense::Vec< T, 3 >  g_lat[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints],
const T  r_inv,
const T  grad_r_inv,
const int  wedge_idx,
const int  node_idx,
const int  quad_point_idx 
)
constexpr

Computes and returns J^-T grad(N_j).

Parameters
g_rad[in] see lateral_parts_of_grad_phi()
g_lat[in] see lateral_parts_of_grad_phi()
r_inv[in] 1 / r (where r is the physical radius of the quadrature point - compute with forward_map_rad())
grad_r_inv[in] 1 / grad_r (where grad_r is the radial component of the gradient of the forward map in radial direction - compute with grad_forward_map_rad())
wedge_idx[in] wedge index of the hex cell (0 or 1)
node_idx[in] node index in the wedge (0, 1, ..., 5)
quad_point_idx[in] index of the quadrature point

◆ grad_shape_lat_coarse_eta()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_shape_lat_coarse_eta ( const int  coarse_node_idx,
const int  fine_lateral_wedge_idx 
)
constexpr

◆ grad_shape_lat_coarse_xi()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_shape_lat_coarse_xi ( const int  coarse_node_idx,
const int  fine_lateral_wedge_idx 
)
constexpr

◆ grad_shape_lat_eta()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_shape_lat_eta ( const int  node_idx)
constexpr

Gradient of the lateral part of the shape function, in eta direction \( \frac{\partial}{\partial \eta} N^\mathrm{lat}_j \).

This is different from the \( \eta \) part (second entry) of the gradient of the full shape function!

That would be

\( \frac{\partial}{\partial \eta} N_j = N^\mathrm{rad}_j \frac{\partial}{\partial \eta} N^\mathrm{lat}_j \)

◆ grad_shape_lat_xi()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_shape_lat_xi ( const int  node_idx)
constexpr

Gradient of the lateral part of the shape function, in xi direction \( \frac{\partial}{\partial \xi} N^\mathrm{lat}_j \).

This is different from the \( \xi \) part (first entry) of the gradient of the full shape function!

That would be

\( \frac{\partial}{\partial \xi} N_j = N^\mathrm{rad}_j \frac{\partial}{\partial \xi} N^\mathrm{lat}_j \)

◆ grad_shape_rad()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_shape_rad ( const int  node_idx)
constexpr

Gradient of the radial part of the shape function, in the radial direction \( \frac{\partial}{\partial \zeta} N^\mathrm{rad}_j \).

This is different from the radial part of the gradient of the full shape function!

That would be

\( \frac{\partial}{\partial \zeta} N_j = N^\mathrm{lat}_j \frac{\partial}{\partial \zeta} N^\mathrm{rad}_j \)

◆ grad_shape_rad_coarse()

template<std::floating_point T>
constexpr T terra::fe::wedge::grad_shape_rad_coarse ( const int  coarse_node_idx,
const int  fine_radial_wedge_idx 
)
constexpr

◆ jac() [1/2]

template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > terra::fe::wedge::jac ( const dense::Vec< T, 3 > &  p1_phy,
const dense::Vec< T, 3 > &  p2_phy,
const dense::Vec< T, 3 > &  p3_phy,
const T  r_1,
const T  r_2,
const T  xi,
const T  eta,
const T  zeta 
)
constexpr

◆ jac() [2/2]

template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > terra::fe::wedge::jac ( const dense::Vec< T, 3 >  p_phy[3],
const T  r_1,
const T  r_2,
const dense::Vec< T, 3 > &  xi_eta_zeta_fine 
)
constexpr

◆ jac_lat()

template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > terra::fe::wedge::jac_lat ( const dense::Vec< T, 3 > &  p1_phy,
const dense::Vec< T, 3 > &  p2_phy,
const dense::Vec< T, 3 > &  p3_phy,
const T  xi,
const T  eta 
)
constexpr

◆ jac_rad()

template<std::floating_point T>
constexpr dense::Vec< T, 3 > terra::fe::wedge::jac_rad ( const T  r_1,
const T  r_2,
const T  zeta 
)
constexpr

◆ jacobian_lat_determinant()

template<std::floating_point T, int NumQuadPoints>
constexpr void terra::fe::wedge::jacobian_lat_determinant ( T(&)  det_jac_lat[num_wedges_per_hex_cell][NumQuadPoints],
const dense::Vec< T, 3 >  wedge_surf_phy_coords[num_wedges_per_hex_cell][num_nodes_per_wedge_surface],
const dense::Vec< T, 3 >  quad_points[NumQuadPoints] 
)
constexpr

Like jacobian_lat_inverse_transposed_and_determinant() but only computes the determinant (cheaper if the transposed inverse of the Jacobian is not required).

Parameters
det_jac_lat[out] absolute of the determinant of the Jacobian of the lateral map
wedge_surf_phy_coords[in] coords of the triangular wedge surfaces on the unit sphere (compute via wedge_surface_physical_coords())
quad_points[in] the quadrature points

◆ jacobian_lat_inverse_transposed_and_determinant()

template<std::floating_point T, int NumQuadPoints>
constexpr void terra::fe::wedge::jacobian_lat_inverse_transposed_and_determinant ( dense::Mat< T, 3, 3 >(&)  jac_lat_inv_t[num_wedges_per_hex_cell][NumQuadPoints],
T(&)  det_jac_lat[num_wedges_per_hex_cell][NumQuadPoints],
const dense::Vec< T, 3 >  wedge_surf_phy_coords[num_wedges_per_hex_cell][num_nodes_per_wedge_surface],
const dense::Vec< T, 3 >  quad_points[NumQuadPoints] 
)
constexpr

Computes the transposed inverse of the Jacobian of the lateral forward map from the reference triangle to the triangle on the unit sphere and the absolute determinant of that Jacobian at the passed quadrature points.

Parameters
jac_lat_inv_t[out] transposed inverse of the Jacobian of the lateral map
det_jac_lat[out] absolute of the determinant of the Jacobian of the lateral map
wedge_surf_phy_coords[in] coords of the triangular wedge surfaces on the unit sphere (compute via wedge_surface_physical_coords())
quad_points[in] the quadrature points

◆ lateral_parts_of_grad_phi()

template<std::floating_point T, int NumQuadPoints>
constexpr void terra::fe::wedge::lateral_parts_of_grad_phi ( dense::Vec< T, 3 >(&)  g_rad[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints],
dense::Vec< T, 3 >(&)  g_lat[num_wedges_per_hex_cell][num_nodes_per_wedge][NumQuadPoints],
const dense::Mat< T, 3, 3 >  jac_lat_inv_t[num_wedges_per_hex_cell][NumQuadPoints],
const dense::Vec< T, 3 >  quad_points[NumQuadPoints] 
)
constexpr

Computes the radially independent parts of the physical shape function gradients.

g_rad_j = jac_lat_inv_t * ( (∂/∂xi N_lat_j) N_rad_j, (∂/∂eta N_lat_j) N_rad_j, 0 )^T
g_lat_j = jac_lat_inv_t * ( 0, 0, N_lat_j )^T

where j is a node of the wedge. Computes those for all 6 nodes j = 0, ..., 5 of a wedge.

Those terms can technically be precomputed and are independent of r (identical for all wedges on a radial beam).

From thes we can later compute

jac_inv_t * grad_N_j

via

jac_inv_t * grad_N_j = (1 / r(zeta)) g_rad_j + (∂ N_rad_j / ∂zeta) * (1 / (∂r / ∂zeta)) * g_lat_j
Parameters
g_rad[out] g_rad - see above
g_lat[out] g_lat - see above
jac_lat_inv_t[in] transposed inverse of lateral Jacobian - see jacobian_lat_inverse_transposed_and_determinant()
quad_points[in] the quadrature points

◆ reorder_local_dofs()

template<typename ScalarT >
constexpr void terra::fe::wedge::reorder_local_dofs ( const DoFOrdering  doo_from,
const DoFOrdering  doo_to,
dense::Vec< ScalarT, 18 > &  dofs 
)
constexpr

◆ shape() [1/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape ( const int  node_idx,
const dense::Vec< T, 3 > &  xi_eta_zeta 
)
constexpr

(Tensor-product) Shape function.

\[ N_i = N^\mathrm{lat}_i * N^\mathrm{rad}_i \]

◆ shape() [2/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape ( const int  node_idx,
const T  xi,
const T  eta,
const T  zeta 
)
constexpr

(Tensor-product) Shape function.

\[ N_i = N^\mathrm{lat}_i * N^\mathrm{rad}_i \]

◆ shape_coarse() [1/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_coarse ( const int  coarse_node_idx,
const int  fine_radial_wedge_idx,
const int  fine_lateral_wedge_idx,
const dense::Vec< T, 3 > &  xi_eta_zeta_fine 
)
constexpr

◆ shape_coarse() [2/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_coarse ( const int  coarse_node_idx,
const int  fine_radial_wedge_idx,
const int  fine_lateral_wedge_idx,
const T  xi_fine,
const T  eta_fine,
const T  zeta_fine 
)
constexpr

◆ shape_lat() [1/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_lat ( const int  node_idx,
const dense::Vec< T, 3 > &  xi_eta_zeta 
)
constexpr

Lateral shape function.

\[ \begin{align} N^\mathrm{lat}_0 = N^\mathrm{lat}_3 &= 1 - \xi - \eta \\ N^\mathrm{lat}_1 = N^\mathrm{lat}_4 &= \xi \\ N^\mathrm{lat}_2 = N^\mathrm{lat}_5 &= \eta \end{align} \]

◆ shape_lat() [2/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_lat ( const int  node_idx,
const T  xi,
const T  eta 
)
constexpr

Lateral shape function.

\[ \begin{align} N^\mathrm{lat}_0 = N^\mathrm{lat}_3 &= 1 - \xi - \eta \\ N^\mathrm{lat}_1 = N^\mathrm{lat}_4 &= \xi \\ N^\mathrm{lat}_2 = N^\mathrm{lat}_5 &= \eta \end{align} \]

◆ shape_lat_coarse()

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_lat_coarse ( const int  coarse_node_idx,
const int  fine_lateral_wedge_idx,
const T  xi_fine,
const T  eta_fine 
)
constexpr

Returns the coarse grid lateral shape function evaluated at a point of the reference fine grid wedge.

Parameters
coarse_node_idxwedge node index of the coarse grid shape function ("which coarse grid shape function to evaluate") (in {0, ..., 5})
fine_lateral_wedge_idxindex of the lateral wedge index in a (once) refined coarse mesh, in {0, 1, 2, 3} 0: bottom left triangle (orientation up) 1: bottom right triangle (orientation: up) 2: top triangle (orientation: up) 3: center triangle (orientation: down)
xi_finexi-coordinate in the reference fine-wedge (in [0, 1])
eta_fineeta-coordinate in the reference fine-wedge (in [0, 1])

◆ shape_rad() [1/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_rad ( const int  node_idx,
const dense::Vec< T, 3 > &  xi_eta_zeta 
)
constexpr

Radial shape function.

\[ \begin{align} N^\mathrm{rad}_0 = N^\mathrm{rad}_1 = N^\mathrm{rad}_2 &= \frac{1}{2} ( 1 - \zeta ) \\ N^\mathrm{rad}_3 = N^\mathrm{rad}_4 = N^\mathrm{rad}_5 &= \frac{1}{2} ( 1 + \zeta ) \\ \end{align} \]

◆ shape_rad() [2/2]

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_rad ( const int  node_idx,
const T  zeta 
)
constexpr

Radial shape function.

\[ \begin{align} N^\mathrm{rad}_0 = N^\mathrm{rad}_1 = N^\mathrm{rad}_2 &= \frac{1}{2} ( 1 - \zeta ) \\ N^\mathrm{rad}_3 = N^\mathrm{rad}_4 = N^\mathrm{rad}_5 &= \frac{1}{2} ( 1 + \zeta ) \\ \end{align} \]

◆ shape_rad_coarse()

template<std::floating_point T>
constexpr T terra::fe::wedge::shape_rad_coarse ( const int  coarse_node_idx,
const int  fine_radial_wedge_idx,
const T  zeta_fine 
)
constexpr

Returns the coarse grid radial shape function evaluated at a point of the reference fine grid wedge.

Parameters
coarse_node_idxwedge node index of the coarse grid shape function ("which coarse grid shape function to evaluate") (in {0, ..., 5})
fine_radial_wedge_idx0 for inner fine wedge, 1 for outer fine wedge
zeta_finecoordinate in the reference fine-wedge (in [-1, 1])

◆ symmetric_grad()

template<std::floating_point T>
constexpr dense::Mat< T, 3, 3 > terra::fe::wedge::symmetric_grad ( const dense::Mat< T, 3, 3 > &  J_inv_transposed,
const dense::Vec< T, 3 > &  quad_point,
const int  dof,
const int  dim 
)
constexpr

Returns the symmetric gradient of the shape function of a dof at a quadrature point.

Parameters
J_inv_transposedinverse transposed jacobian to map to physical element
quad_pointquadrature point on the reference element
doflocal index of the shape function a wedge
dimdimension of the vectorial shape function, of which the gradient should be computed

◆ wedge_0_surface_physical_coords()

template<std::floating_point T>
void terra::fe::wedge::wedge_0_surface_physical_coords ( dense::Vec< T, 3 > *  wedge_surf_phy_coords,
const grid::Grid3DDataVec< T, 3 > &  lateral_grid,
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell 
)

◆ wedge_1_surface_physical_coords()

template<std::floating_point T>
void terra::fe::wedge::wedge_1_surface_physical_coords ( dense::Vec< T, 3 > *  wedge_surf_phy_coords,
const grid::Grid3DDataVec< T, 3 > &  lateral_grid,
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell 
)

◆ wedge_surface_physical_coords()

template<std::floating_point T>
void terra::fe::wedge::wedge_surface_physical_coords ( dense::Vec< T, 3 >(&)  wedge_surf_phy_coords[num_wedges_per_hex_cell][num_nodes_per_wedge_surface],
const grid::Grid3DDataVec< T, 3 > &  lateral_grid,
const int  local_subdomain_id,
const int  x_cell,
const int  y_cell 
)

Extracts the (unit sphere) surface vertex coords of the two wedges of a hex cell.

Useful for wedge-based kernels that update two wedges that make up a hex cell at once.

2--3
|\ |
| \| => [(p0, p1, p2), (p3, p2, p1)]
0--1
Parameters
wedge_surf_phy_coords[out] first dim: wedge/triangle index, second dim: vertex index
lateral_grid[in] the unit sphere vertex coordinates
local_subdomain_id[in] shell subdomain id on this process
x_cell[in] hex cell x-coordinate
y_cell[in] hex cell y-coordinate

Variable Documentation

◆ num_nodes_per_wedge

constexpr int terra::fe::wedge::num_nodes_per_wedge = 6
constexpr

◆ num_nodes_per_wedge_surface

constexpr int terra::fe::wedge::num_nodes_per_wedge_surface = 3
constexpr

◆ num_wedges_per_hex_cell

constexpr int terra::fe::wedge::num_wedges_per_hex_cell = 2
constexpr