NodesAndModes
NodesAndModes.abctorst
— Methodabctorst(elem::Pyr,a,b,c)
Converts from Stroud coordinates (a,b,c) on [-1,1]^3 to reference element coordinates (r,s,t).
NodesAndModes.abtors
— Methodabtors(Tri(), r, s, tol = 1e-12)
Converts from polynomial basis evaluation coordinates (a,b) on the domain [-1,1]^2 to reference bi-unit right triangle coordinate (r,s).
NodesAndModes.basis
— Functionbasis(elem::Pyr,N,r,s,t,tol=1e-12)
Computes orthonormal semi-nodal basis on the biunit pyramid element.
Warning: nodal derivative matrices may contain errors for nodes at t = 1. A way to avoid this is to use weak differentiation matrices computed using quadrature rules with only interior nodes.
NodesAndModes.basis
— Methodbasis(elem::AbstractElemShape, N, rst...)
Computes orthonormal basis of degree N at tuple of coordinate arrays (r,s,t).
NodesAndModes.basis
— Methodbasis(elem::Line,N,r)
Computes the generalized Vandermonde matrix V of degree N (along with the derivative matrix Vr) at points r.
NodesAndModes.build_warped_nodes
— Methodbuild_warped_nodes(elem::AbstractElemShape, N, r1D)
Computes degree N warp-and-blend interpolation nodes for elem = Tri(), Pyr(), or Tet() based on the 1D node set "r1D". Returns a tuple "rst" containing arrays of interpolation points.
NodesAndModes.edge_basis
— Methodedge_basis(elem::AbstractElemShape, N, rst...)
Returns the generalized Vandermonde matrix evaluated using an edge basis (e.g., degree N
polynomials over an edge, but linearly blended into the interior). The dimension of the resulting space is simply the number of total nodes on edges of a degree N
element.
NodesAndModes.edge_basis
— Methodedge_basis(N, vertices, edges, basis, vertex_functions, rst...)
Computes edge basis given vertex functions and 1D basis.
NodesAndModes.equi_nodes
— Methodequi_nodes(elem::AbstractElemShape, N)
Compute equispaced nodes of degree N.
NodesAndModes.equi_nodes
— Methodequi_nodes(elem::Line,N)
Computes equally spaced nodes of degree N.
NodesAndModes.find_face_nodes
— Functionfunction find_face_nodes(elem, r, s, tol=50*eps())
function find_face_nodes(elem, r, s, t, tol=50*eps())
Given volume nodes r
, s
, t
, finds face nodes. Note that this function implicitly defines an ordering on the faces.
NodesAndModes.gauss_lobatto_quad
— Methodgauss_lobatto_quad(α, β, N)
Computes Legendre-Gauss-Lobatto quadrature points and weights with Jacobi weights α,β.
NodesAndModes.gauss_quad
— Methodgauss_quad(α, β, N)
Compute nodes and weights for Gaussian quadrature with Jacobi weights (α,β).
NodesAndModes.get_edge_list
— Methodget_edge_list(elem::AbstractElemShape)
Returns list of edges for a specific element (elem = Tri(), Pyr(), Hex(), or Tet()).
NodesAndModes.grad_jacobiP
— Methodgrad_jacobiP(r, α, β, N, tmp_array=())
Evaluate derivative of Jacobi polynomial (α, β) of order N at points r.
NodesAndModes.grad_vandermonde
— Methodgrad_vandermonde(elem::AbstractElemShape, N, rst...)
Computes the generalized Vandermonde derivative matrix V of degree N at points (r,s,t).
NodesAndModes.interp_1D_to_edges
— Methodinterp_1D_to_edges(elem::AbstractElemShape, r1D)
Interpolates points r1D to the edges of an element (elem = :Tri, :Pyr, or :Tet)
NodesAndModes.jaskowiec_sukumar_quad_nodes
— Methodjaskowiek_sukumar_quad_nodes(elem::Tet)
Symmetric quadrature rules on the tetrahedron of degree up to 20 from:
Jaśkowiec, J, Sukumar, N., "High-order symmetric cubature rules for tetrahedra and pyramids." Int J Numer Methods Eng. 122(1): 148-171, 2021.
NodesAndModes.meshgrid
— Methodmeshgrid(vx) Computes an (x,y)-grid from the vectors (vx,vx). For more information, see the MATLAB documentation.
Copied and pasted directly from VectorizedRoutines.jl. Using VectorizedRoutines.jl directly causes Pkg versioning issues with SpecialFunctions.jl
NodesAndModes.meshgrid
— Methodmeshgrid(vx,vy,vz) Computes an (x,y,z)-grid from the vectors (vx,vy,vz). For more information, see the MATLAB documentation.
Copied and pasted directly from VectorizedRoutines.jl. Using VectorizedRoutines.jl directly causes Pkg versioning issues with SpecialFunctions.jl
NodesAndModes.meshgrid
— Methodmeshgrid(vx,vy) Computes an (x,y)-grid from the vectors (vx,vy). For more information, see the MATLAB documentation.
Copied and pasted directly from VectorizedRoutines.jl. Using VectorizedRoutines.jl directly causes Pkg versioning issues with SpecialFunctions.jl
NodesAndModes.nodes
— Methodnodes(elem::AbstractElemShape,N)
Computes interpolation nodes of degree N. Edge nodes coincide with (N+1)-point Lobatto points. Default routine for elem = Tet(), Pyr(), Tri().
For Quad(), Hex(), Wedge() elements, nodes(...) returns interpolation points constructed using a tensor product of lower-dimensional nodes.
NodesAndModes.nodes
— Methodnodes(elem::Line,N)
Computes interpolation nodes of degree N.
NodesAndModes.nodes
— Methodnodes(elem::Pyr,N)
Computes interpolation nodes of degree N. Edge nodes coincide with (N+1)-point Lobatto points. Triangular face nodes coincide with Tri.nodes(N), quadrilateral face nodes coincide with tensor product (N+1)-point Lobatto points.
NodesAndModes.quad_nodes
— Methodquad_nodes(elem::AbstractElemShape, N)
Compute quadrature nodes and weights exact for (at least) degree 2N polynomials.
NodesAndModes.quad_nodes
— Methodquad_nodes(elem::Line,N)
Computes (N+1)-point Gauss quadrature rule (exact for degree 2N+1 polynomials)
NodesAndModes.quad_nodes_tet
— Methodquad_nodes_tet(N)
Returns quadrature nodes and weights which exactly integrate degree N polynomials
NodesAndModes.quad_nodes_tri
— Methodquad_nodes_tri(N)
Returns quadrature nodes (from Gimbutas and Xiao 2010) which exactly integrate degree N polynomials
NodesAndModes.rstoab
— Functionrstoab(r, s, tol = 1e-12)
Converts from reference bi-unit right triangle coordinate (r,s) to polynomial basis evaluation coordinates (a,b) on the domain [-1,1]^2
NodesAndModes.rsttoabc
— Methodrsttoabc(elem::Pyr,a,b,c)
Converts from reference element coordinates (r,s,t) to Stroud coordinates (a,b,c) on [-1,1]^3.
NodesAndModes.simplex_2D
— Methodsimplex_2D(a, b, i, j)
Evaluate 2D PKDO basis phi_ij at points (a,b) on the Duffy domain [-1,1]^2
NodesAndModes.simplex_3D
— Methodsimplex_3D(a, b, c, i, j, k)
Evaluate 3D "Legendre" basis phi_ijk at (a,b,c) coordinates on the [-1,1] cube
NodesAndModes.stroud_quad_nodes
— Methodstroud_quad_nodes(elem::AbstractElemShape,N)
Returns Stroud-type quadrature nodes and weights constructed from the tensor product of (N+1)-point Gauss-Jacobi rules. Exact for degree 2N polynomials
NodesAndModes.vandermonde
— Methodvandermonde(elem::AbstractElemShape, N, rst...)
Computes the generalized Vandermonde matrix V of degree N at points (r,s,t).