Go to the documentation of this file.
95 #ifndef INTREPID_POLYLIB_HPP
96 #define INTREPID_POLYLIB_HPP
98 #include "Intrepid_ConfigDefs.hpp"
100 #include "Teuchos_Assert.hpp"
200 PL_GAUSS_RADAU_RIGHT,
205 inline EIntrepidPLPoly & operator++(EIntrepidPLPoly &type) {
206 return type = static_cast<EIntrepidPLPoly>(type+1);
236 template<
class Scalar>
237 static void zwgj (Scalar *z, Scalar *w,
const int np,
const Scalar alpha,
const Scalar beta);
247 template<
class Scalar>
248 static void zwgrjm (Scalar *z, Scalar *w,
const int np,
const Scalar alpha,
const Scalar beta);
258 template<
class Scalar>
259 static void zwgrjp (Scalar *z, Scalar *w,
const int np,
const Scalar alpha,
const Scalar beta);
269 template<
class Scalar>
270 static void zwglj (Scalar *z, Scalar *w,
const int np,
const Scalar alpha,
const Scalar beta);
284 template<
class Scalar>
285 static void Dgj (Scalar *D,
const Scalar *z,
const int np,
const Scalar alpha,
const Scalar beta);
296 template<
class Scalar>
297 static void Dgrjm (Scalar *D,
const Scalar *z,
const int np,
const Scalar alpha,
const Scalar beta);
308 template<
class Scalar>
309 static void Dgrjp (Scalar *D,
const Scalar *z,
const int np,
const Scalar alpha,
const Scalar beta);
320 template<
class Scalar>
321 static void Dglj (Scalar *D,
const Scalar *z,
const int np,
const Scalar alpha,
const Scalar beta);
346 template<
class Scalar>
347 static Scalar
hgj (
const int i,
const Scalar z,
const Scalar *zgj,
348 const int np,
const Scalar alpha,
const Scalar beta);
370 template<
class Scalar>
371 static Scalar
hgrjm (
const int i,
const Scalar z,
const Scalar *zgrj,
372 const int np,
const Scalar alpha,
const Scalar beta);
394 template<
class Scalar>
395 static Scalar
hgrjp (
const int i,
const Scalar z,
const Scalar *zgrj,
396 const int np,
const Scalar alpha,
const Scalar beta);
418 template<
class Scalar>
419 static Scalar
hglj (
const int i,
const Scalar z,
const Scalar *zglj,
420 const int np,
const Scalar alpha,
const Scalar beta);
436 template<
class Scalar>
437 static void Imgj (Scalar *im,
const Scalar *zgj,
const Scalar *zm,
const int nz,
438 const int mz,
const Scalar alpha,
const Scalar beta);
451 template<
class Scalar>
452 static void Imgrjm(Scalar *im,
const Scalar *zgrj,
const Scalar *zm,
const int nz,
453 const int mz,
const Scalar alpha,
const Scalar beta);
466 template<
class Scalar>
467 static void Imgrjp(Scalar *im,
const Scalar *zgrj,
const Scalar *zm,
const int nz,
468 const int mz,
const Scalar alpha,
const Scalar beta);
481 template<
class Scalar>
482 static void Imglj (Scalar *im,
const Scalar *zglj,
const Scalar *zm,
const int nz,
483 const int mz,
const Scalar alpha,
const Scalar beta);
527 template<
class Scalar>
528 static void jacobfd (
const int np,
const Scalar *z, Scalar *poly_in, Scalar *polyd,
529 const int n,
const Scalar alpha,
const Scalar beta);
545 template<
class Scalar>
546 static void jacobd (
const int np,
const Scalar *z, Scalar *polyd,
const int n,
547 const Scalar alpha,
const Scalar beta);
559 template<
class Scalar>
560 static void Jacobz (
const int n, Scalar *z,
const Scalar alpha,
const Scalar beta);
585 template<
class Scalar>
586 static void JacZeros (
const int n, Scalar *a,
const Scalar alpha,
const Scalar beta);
612 template<
class Scalar>
613 static void TriQL (
const int n, Scalar *d, Scalar *e);
625 template<
class Scalar>
626 static Scalar
gammaF (
const Scalar x);
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.
static void zwgj(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Jacobi zeros and weights.
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...
Contains definitions of custom data types in Intrepid.
EIntrepidPLPoly
Enumeration of coordinate frames (reference/ambient) for geometrical entities (cells,...
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...
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.
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.
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.
Providing orthogonal polynomial calculus and interpolation, created by Spencer Sherwin,...
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.
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.
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 ...
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 ar...
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.
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. .
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...
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.
static void TriQL(const int n, Scalar *d, Scalar *e)
QL algorithm for symmetric tridiagonal matrix.
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 rela...
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 ...
static Scalar gammaF(const Scalar x)
Calculate the Gamma function , , for integer values x and halves.
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 z...
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.
Definition file for a set of functions providing orthogonal polynomial polynomial calculus and interp...
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, .