|
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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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, . | |
| 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. | |
| 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. . | |
| 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. | |
| template<class Scalar > | |
| static void | TriQL (const int n, Scalar *d, Scalar *e) |
| QL algorithm for symmetric tridiagonal matrix. | |
| template<class Scalar > | |
| static Scalar | gammaF (const Scalar x) |
Calculate the Gamma function , , for integer values x and halves. | |
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.
| void Intrepid::IntrepidPolylib::Dgj | ( | Scalar * | D, |
| const Scalar * | z, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::Dglj | ( | Scalar * | D, |
| const Scalar * | z, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [static] |
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros.
Definition at line 328 of file Intrepid_PolylibDef.hpp.
| void Intrepid::IntrepidPolylib::Dgrjm | ( | Scalar * | D, |
| const Scalar * | z, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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.
| void Intrepid::IntrepidPolylib::Dgrjp | ( | Scalar * | D, |
| const Scalar * | z, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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.
| Scalar Intrepid::IntrepidPolylib::gammaF | ( | const Scalar | x | ) | [static] |
Calculate the Gamma function ,
, for integer values x and halves.
Determine the value of
using:

where
Definition at line 806 of file Intrepid_PolylibDef.hpp.
Referenced by Dglj(), Dgrjm(), Dgrjp(), JacZeros(), zwgj(), zwglj(), zwgrjm(), and zwgrjp().
| Scalar Intrepid::IntrepidPolylib::hgj | ( | const int | i, |
| const Scalar | z, | ||
| const Scalar * | zgj, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| Scalar Intrepid::IntrepidPolylib::hglj | ( | const int | i, |
| const Scalar | z, | ||
| const Scalar * | zglj, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| Scalar Intrepid::IntrepidPolylib::hgrjm | ( | const int | i, |
| const Scalar | z, | ||
| const Scalar * | zgrj, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| Scalar Intrepid::IntrepidPolylib::hgrjp | ( | const int | i, |
| const Scalar | z, | ||
| const Scalar * | zgrj, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::Imgj | ( | Scalar * | im, |
| const Scalar * | zgj, | ||
| const Scalar * | zm, | ||
| const int | nz, | ||
| const int | mz, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::Imglj | ( | Scalar * | im, |
| const Scalar * | zglj, | ||
| const Scalar * | zm, | ||
| const int | nz, | ||
| const int | mz, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::Imgrjm | ( | Scalar * | im, |
| const Scalar * | zgrj, | ||
| const Scalar * | zm, | ||
| const int | nz, | ||
| const int | mz, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::Imgrjp | ( | Scalar * | im, |
| const Scalar * | zgrj, | ||
| const Scalar * | zm, | ||
| const int | nz, | ||
| const int | mz, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::jacobd | ( | const int | np, |
| const Scalar * | z, | ||
| Scalar * | polyd, | ||
| const int | n, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [static] |
Calculate the derivative of Jacobi polynomials.
at the np points z.
Definition at line 653 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().
| void Intrepid::IntrepidPolylib::jacobfd | ( | const int | np, |
| const Scalar * | z, | ||
| Scalar * | poly_in, | ||
| Scalar * | polyd, | ||
| const int | n, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [static] |
Routine to calculate Jacobi polynomials,
, and their first derivative,
.
and its derivative at the np points in z[i]![$ \begin{array}{rcl} P^{\alpha,\beta}_0(z) &=& 1 \\ P^{\alpha,\beta}_1(z) &=& \frac{1}{2} [ \alpha-\beta+(\alpha+\beta+2)z] \\ a^1_n P^{\alpha,\beta}_{n+1}(z) &=& (a^2_n + a^3_n z) P^{\alpha,\beta}_n(z) - a^4_n P^{\alpha,\beta}_{n-1}(z) \\ a^1_n &=& 2(n+1)(n+\alpha + \beta + 1)(2n + \alpha + \beta) \\ a^2_n &=& (2n + \alpha + \beta + 1)(\alpha^2 - \beta^2) \\ a^3_n &=& (2n + \alpha + \beta)(2n + \alpha + \beta + 1) (2n + \alpha + \beta + 2) \\ a^4_n &=& 2(n+\alpha)(n+\beta)(2n + \alpha + \beta + 2) \end{array} $](form_138.png)
![$ \begin{array}{rcl} b^1_n(z)\frac{d}{dz} P^{\alpha,\beta}_n(z)&=&b^2_n(z)P^{\alpha,\beta}_n(z) + b^3_n(z) P^{\alpha,\beta}_{n-1}(z) \hspace{2.2cm} \\ b^1_n(z) &=& (2n+\alpha + \beta)(1-z^2) \\ b^2_n(z) &=& n[\alpha - \beta - (2n+\alpha + \beta)z]\\ b^3_n(z) &=& 2(n+\alpha)(n+\beta) \end{array} $](form_139.png)
Definition at line 568 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().
| void Intrepid::IntrepidPolylib::Jacobz | ( | const int | n, |
| Scalar * | z, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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 670 of file Intrepid_PolylibDef.hpp.
References INTREPID_POLYLIB_STOP, and jacobfd().
| void Intrepid::IntrepidPolylib::JacZeros | ( | const int | n, |
| Scalar * | a, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [static] |
Zero determination through the eigenvalues of a tridiagonal matrix from the three term recursion relationship.
Set up a symmetric tridiagonal matrix
![$ \left [ \begin{array}{ccccc} a[0] & b[0] & & & \\ b[0] & a[1] & b[1] & & \\ 0 & \ddots & \ddots & \ddots & \\ & & \ddots & \ddots & b[n-2] \\ & & & b[n-2] & a[n-1] \end{array} \right ] $](form_144.png)
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 701 of file Intrepid_PolylibDef.hpp.
| void Intrepid::IntrepidPolylib::TriQL | ( | const int | n, |
| Scalar * | d, | ||
| Scalar * | e | ||
| ) | [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 740 of file Intrepid_PolylibDef.hpp.
Referenced by JacZeros().
| void Intrepid::IntrepidPolylib::zwgj | ( | Scalar * | z, |
| Scalar * | w, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::zwglj | ( | Scalar * | z, |
| Scalar * | w, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::zwgrjm | ( | Scalar * | z, |
| Scalar * | w, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
| void Intrepid::IntrepidPolylib::zwgrjp | ( | Scalar * | z, |
| Scalar * | w, | ||
| const int | np, | ||
| const Scalar | alpha, | ||
| const Scalar | beta | ||
| ) | [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().
1.7.6.1