Defines the reference tetrahedron element used in nodal Discontinuous Galerkin methods. More...
#include <xtensor/core/xtensor_forward.hpp>#include "oiseau/dg/nodal/ref_element.hpp"Go to the source code of this file.
Classes | |
| class | oiseau::dg::nodal::RefTetrahedron |
| Represents a reference tetrahedron element for nodal Discontinuous Galerkin methods. More... | |
Functions | |
| xt::xarray< double > | oiseau::dg::nodal::detail::equilateral_xyz_to_rst (const xt::xarray< double > &coords) |
| Transforms coordinates from a canonical "equilateral" tetrahedron in XYZ space to the reference (r,s,t) computational space. | |
| xt::xarray< double > | oiseau::dg::nodal::detail::rst_to_abc (const xt::xarray< double > &rst) |
| Transforms coordinates from the reference (r,s,t) space to a auxiliary (a,b,c) space. | |
| xt::xarray< double > | oiseau::dg::nodal::detail::generate_tetrahedron_nodes (unsigned p) |
| Generates nodal points on the reference tetrahedron for a given polynomial order. | |
| xt::xarray< double > | oiseau::dg::nodal::detail::generate_tetrahedron_equidistant_nodes (unsigned n) |
| Generates approximately equidistant nodes on the reference tetrahedron. | |
| xt::xarray< double > | oiseau::dg::nodal::detail::warp_shift_face_3d (int p, double pval, const xt::xarray< double > &l2, const xt::xarray< double > &l3, const xt::xarray< double > &l4) |
| Helper function for the Warp & Blend algorithm: computes node shifts for points near the faces of the tetrahedron. | |
| xt::xarray< double > | oiseau::dg::nodal::detail::eval_shift (int p, double pval, const xt::xarray< double > &l1, const xt::xarray< double > &l2, const xt::xarray< double > &l3) |
| Helper function for the Warp & Blend algorithm: evaluates a shifting polynomial. | |
| xt::xarray< double > | oiseau::dg::nodal::detail::eval_warp (int p, const xt::xarray< double > &xnodes, const xt::xarray< double > &xout) |
| Helper function for the Warp & Blend algorithm: evaluates the warp factor or displacement. | |
Defines the reference tetrahedron element used in nodal Discontinuous Galerkin methods.
This file contains the declaration of the RefTetrahedron class, which represents a reference tetrahedron element, and associated internal helper functions in the detail namespace for generating nodes and performing coordinate transformations.
| xt::xarray< double > oiseau::dg::nodal::detail::equilateral_xyz_to_rst | ( | const xt::xarray< double > & | coords | ) |
Transforms coordinates from a canonical "equilateral" tetrahedron in XYZ space to the reference (r,s,t) computational space.
The "equilateral" tetrahedron might be one with vertices like (0,0,0), (1,0,0), (0.5, sqrt(3)/2, 0), (0.5, sqrt(3)/6, sqrt(6)/3). The reference (r,s,t) space is typically { (r,s,t) | r,s,t >= -1, r+s+t <= -1 }.
| coords | 2D array (shape: N_points × 3) of XYZ coordinates. |
| xt::xarray< double > oiseau::dg::nodal::detail::eval_shift | ( | int | p, |
| double | pval, | ||
| const xt::xarray< double > & | l1, | ||
| const xt::xarray< double > & | l2, | ||
| const xt::xarray< double > & | l3 ) |
Helper function for the Warp & Blend algorithm: evaluates a shifting polynomial.
This polynomial contributes to the overall warp transformation, often based on distances to element boundaries (represented by barycentric coordinates).
| p | Polynomial order. |
| pval | A parameter, possibly related to the evaluation point or a scaling factor for the shift. |
| l1 | Barycentric coordinate \(\zeta_1\). |
| l2 | Barycentric coordinate \(\zeta_2\). |
| l3 | Barycentric coordinate \(\zeta_3\). |
| xt::xarray< double > oiseau::dg::nodal::detail::eval_warp | ( | int | p, |
| const xt::xarray< double > & | xnodes, | ||
| const xt::xarray< double > & | xout ) |
Helper function for the Warp & Blend algorithm: evaluates the warp factor or displacement.
This function computes the actual warp displacement for a set of points xout based on the polynomial order p and potentially a reference set of 1D nodes xnodes (e.g., alpha-optimized Gauss-Lobatto nodes).
| p | Polynomial order. |
| xnodes | 1D array of reference nodes (e.g., roots of a Jacobi polynomial) used to define the warp. |
| xout | 1D array of points (e.g., barycentric coordinates) at which to evaluate the warp factor. |
| xt::xarray< double > oiseau::dg::nodal::detail::generate_tetrahedron_equidistant_nodes | ( | unsigned | n | ) |
Generates approximately equidistant nodes on the reference tetrahedron.
These nodes are typically used as a starting point for warping algorithms (like Warp & Blend) or for low-order approximations.
| n | The polynomial order, determining the density of points. The number of points is (n+1)(n+2)(n+3)/6. |
| xt::xarray< double > oiseau::dg::nodal::detail::generate_tetrahedron_nodes | ( | unsigned | p | ) |
Generates nodal points on the reference tetrahedron for a given polynomial order.
These nodes are typically optimized for polynomial interpolation, often generated using methods like Warp & Blend (Warburton, 2006).
| p | The polynomial order. The number of nodes will be (p+1)(p+2)(p+3)/6. |
| xt::xarray< double > oiseau::dg::nodal::detail::rst_to_abc | ( | const xt::xarray< double > & | rst | ) |
Transforms coordinates from the reference (r,s,t) space to a auxiliary (a,b,c) space.
The (a,b,c) coordinates are often collapsed coordinates (e.g., Duffy-like or Koornwinder-type) in which the basis functions (e.g., products of Jacobi polynomials) are conveniently defined. For example, for a tetrahedron, these might be: \(a = r\) \(b = \frac{2(1+s)}{1-t} - 1\) (if t != 1) \(c = t\) (This specific example is for a pyramid, tetrahedral transformations are similar).
| rst | 2D array (shape: N_points × 3) of (r,s,t) coordinates. |
| xt::xarray< double > oiseau::dg::nodal::detail::warp_shift_face_3d | ( | int | p, |
| double | pval, | ||
| const xt::xarray< double > & | l2, | ||
| const xt::xarray< double > & | l3, | ||
| const xt::xarray< double > & | l4 ) |
Helper function for the Warp & Blend algorithm: computes node shifts for points near the faces of the tetrahedron.
This function is part of the process to improve the quality of nodal distributions on the tetrahedron by "warping" an initial set of points.
| p | Polynomial order. |
| pval | A parameter, possibly related to the evaluation point's properties or a scaling factor. |
| l2 | Barycentric coordinate \(\zeta_2\). |
| l3 | Barycentric coordinate \(\zeta_3\). |
| l4 | Barycentric coordinate \(\zeta_4\). |