Intrepid
|
Providing orthogonal polynomial calculus and interpolation, created by Spencer Sherwin, Aeronautics, Imperial College London, modified and redistributed by D. Ridzal. More...
#include <Intrepid_Polylib.hpp>
Static Public Member Functions | |
template<class Scalar > | |
static void | zwgj (Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta) |
Gauss-Jacobi zeros and weights. More... | |
template<class Scalar > | |
static void | zwgrjm (Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta) |
Gauss-Radau-Jacobi zeros and weights with end point at z=-1. More... | |
template<class Scalar > | |
static void | zwgrjp (Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta) |
Gauss-Radau-Jacobi zeros and weights with end point at z=1. More... | |
template<class Scalar > | |
static void | zwglj (Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta) |
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1. More... | |
template<class Scalar > | |
static void | Dgj (Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta) |
Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros. More... | |
template<class Scalar > | |
static void | Dgrjm (Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta) |
Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a zero at z=-1. More... | |
template<class Scalar > | |
static void | Dgrjp (Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta) |
Compute the Derivative Matrix associated with the Gauss-Radau-Jacobi zeros with a zero at z=1. More... | |
template<class Scalar > | |
static void | Dglj (Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta) |
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros. More... | |
template<class Scalar > | |
static Scalar | hgj (const int i, const Scalar z, const Scalar *zgj, const int np, const Scalar alpha, const Scalar beta) |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the arbitrary location z. More... | |
template<class Scalar > | |
static Scalar | hgrjm (const int i, const Scalar z, const Scalar *zgrj, const int np, const Scalar alpha, const Scalar beta) |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point -1. More... | |
template<class Scalar > | |
static Scalar | hgrjp (const int i, const Scalar z, const Scalar *zgrj, const int np, const Scalar alpha, const Scalar beta) |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point +1. More... | |
template<class Scalar > | |
static Scalar | hglj (const int i, const Scalar z, const Scalar *zglj, const int np, const Scalar alpha, const Scalar beta) |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zglj at the arbitrary location z. More... | |
template<class Scalar > | |
static void | Imgj (Scalar *im, const Scalar *zgj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta) |
Interpolation Operator from Gauss-Jacobi points to an arbitrary distribution at points zm. More... | |
template<class Scalar > | |
static void | Imgrjm (Scalar *im, const Scalar *zgrj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta) |
Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at points zm. More... | |
template<class Scalar > | |
static void | Imgrjp (Scalar *im, const Scalar *zgrj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta) |
Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at points zm. More... | |
template<class Scalar > | |
static void | Imglj (Scalar *im, const Scalar *zglj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta) |
Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm. More... | |
template<class Scalar > | |
static void | jacobfd (const int np, const Scalar *z, Scalar *poly_in, Scalar *polyd, const int n, const Scalar alpha, const Scalar beta) |
Routine to calculate Jacobi polynomials, , and their first derivative, . More... | |
template<class Scalar > | |
static void | jacobd (const int np, const Scalar *z, Scalar *polyd, const int n, const Scalar alpha, const Scalar beta) |
Calculate the derivative of Jacobi polynomials. More... | |
template<class Scalar > | |
static void | Jacobz (const int n, Scalar *z, const Scalar alpha, const Scalar beta) |
Calculate the n zeros, z, of the Jacobi polynomial, i.e. . More... | |
template<class Scalar > | |
static void | JacZeros (const int n, Scalar *a, const Scalar alpha, const Scalar beta) |
Zero determination through the eigenvalues of a tridiagonal matrix from the three term recursion relationship. More... | |
template<class Scalar > | |
static void | TriQL (const int n, Scalar *d, Scalar *e) |
QL algorithm for symmetric tridiagonal matrix. More... | |
template<class Scalar > | |
static Scalar | gammaF (const Scalar x) |
Calculate the Gamma function , , for integer values x and halves. More... | |
Providing orthogonal polynomial calculus and interpolation, created by Spencer Sherwin, Aeronautics, Imperial College London, modified and redistributed by D. Ridzal.
See original Polylib documentation.
Definition at line 223 of file Intrepid_Polylib.hpp.
|
static |
Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.
Definition at line 214 of file Intrepid_PolylibDef.hpp.
References jacobd().
|
static |
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros.
Definition at line 328 of file Intrepid_PolylibDef.hpp.
|
static |
Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a zero at z=-1.
Definition at line 247 of file Intrepid_PolylibDef.hpp.
|
static |
Compute the Derivative Matrix associated with the Gauss-Radau-Jacobi zeros with a zero at z=1.
Definition at line 287 of file Intrepid_PolylibDef.hpp.
|
static |
Calculate the Gamma function , , for integer values x and halves.
Determine the value of using:
where
Definition at line 814 of file Intrepid_PolylibDef.hpp.
Referenced by Dglj(), Dgrjm(), Dgrjp(), JacZeros(), zwgj(), zwglj(), zwgrjm(), and zwgrjp().
|
static |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the arbitrary location z.
Definition at line 371 of file Intrepid_PolylibDef.hpp.
References jacobd(), and jacobfd().
Referenced by Imgj().
|
static |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zglj at the arbitrary location z.
Definition at line 434 of file Intrepid_PolylibDef.hpp.
References jacobd(), and jacobfd().
Referenced by Imglj().
|
static |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point -1.
Definition at line 390 of file Intrepid_PolylibDef.hpp.
References jacobd(), and jacobfd().
Referenced by Imgrjm().
|
static |
Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point +1.
Definition at line 412 of file Intrepid_PolylibDef.hpp.
References jacobd(), and jacobfd().
Referenced by Imgrjp().
|
static |
Interpolation Operator from Gauss-Jacobi points to an arbitrary distribution at points zm.
Definition at line 456 of file Intrepid_PolylibDef.hpp.
References hgj().
|
static |
Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm.
Definition at line 540 of file Intrepid_PolylibDef.hpp.
References hglj().
|
static |
Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at points zm.
Definition at line 484 of file Intrepid_PolylibDef.hpp.
References hgrjm().
|
static |
Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at points zm.
Definition at line 512 of file Intrepid_PolylibDef.hpp.
References hgrjp().
|
static |
Calculate the derivative of Jacobi polynomials.
Definition at line 661 of file Intrepid_PolylibDef.hpp.
References jacobfd().
Referenced by Dgj(), Dglj(), Dgrjm(), Dgrjp(), Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >::getValues(), hgj(), hglj(), hgrjm(), hgrjp(), and zwgj().
|
static |
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition at line 570 of file Intrepid_PolylibDef.hpp.
Referenced by Intrepid::Basis_HGRAD_LINE_Cn_FEM_JACOBI< Scalar, ArrayScalar >::getValues(), hgj(), hglj(), hgrjm(), hgrjp(), jacobd(), Jacobz(), zwglj(), zwgrjm(), and zwgrjp().
|
static |
Calculate the n zeros, z, of the Jacobi polynomial, i.e. .
This routine is only valid for and uses polynomial deflation in a Newton iteration
Definition at line 678 of file Intrepid_PolylibDef.hpp.
References INTREPID_POLYLIB_STOP, and jacobfd().
|
static |
Zero determination through the eigenvalues of a tridiagonal matrix from the three term recursion relationship.
Set up a symmetric tridiagonal matrix
Where the coefficients a[n], b[n] come from the recurrence relation
where and are the Jacobi (normalized) orthogonal polynomials ( integer values and halves). Since the polynomials are orthonormalized, the tridiagonal matrix is guaranteed to be symmetric. The eigenvalues of this matrix are the zeros of the Jacobi polynomial.
Definition at line 709 of file Intrepid_PolylibDef.hpp.
|
static |
QL algorithm for symmetric tridiagonal matrix.
This subroutine is a translation of an algol procedure, num. math. 12, 377-383(1968) by martin and wilkinson, as modified in num. math. 15, 450(1970) by dubrulle. Handbook for auto. comp., vol.ii-linear algebra, 241-248(1971). This is a modified version from numerical recipes.
This subroutine finds the eigenvalues and first components of the eigenvectors of a symmetric tridiagonal matrix by the implicit QL method.
on input:
on output:
Definition at line 748 of file Intrepid_PolylibDef.hpp.
Referenced by JacZeros().
|
static |
Gauss-Jacobi zeros and weights.
Definition at line 117 of file Intrepid_PolylibDef.hpp.
References gammaF(), jacobd(), and jacobz.
Referenced by Intrepid::CubaturePolylib< Scalar, ArrayPoint, ArrayWeight >::getCubature(), and Intrepid::PointTools::getGaussPoints().
|
static |
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition at line 186 of file Intrepid_PolylibDef.hpp.
References gammaF(), jacobfd(), and jacobz.
Referenced by Intrepid::CubaturePolylib< Scalar, ArrayPoint, ArrayWeight >::getCubature(), and Intrepid::PointTools::getWarpBlendLatticeLine().
|
static |
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition at line 134 of file Intrepid_PolylibDef.hpp.
References gammaF(), jacobfd(), and jacobz.
Referenced by Intrepid::CubaturePolylib< Scalar, ArrayPoint, ArrayWeight >::getCubature().
|
static |
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition at line 160 of file Intrepid_PolylibDef.hpp.
References gammaF(), jacobfd(), and jacobz.
Referenced by Intrepid::CubaturePolylib< Scalar, ArrayPoint, ArrayWeight >::getCubature().