Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Polylib.h
Go to the documentation of this file.
1 #ifndef H_PLYLIB
2 
5 #include <complex>
6 
7 /*
8  * LIBRARY ROUTINES FOR POLYNOMIAL CALCULUS AND INTERPOLATION
9  */
10 
11 /** \brief The namespace associated with the the Polylib library
12  * (\ref pagePolylib "Polylib introduction")
13  */
14 namespace Polylib {
15 
16  /**
17  \page pagePolylib The Polylib library
18  \section sectionPolyLib Routines For Orthogonal Polynomial Calculus and Interpolation
19 
20  Spencer Sherwin,
21  Aeronautics, Imperial College London
22 
23  Based on codes by Einar Ronquist and Ron Henderson
24 
25  Abbreviations
26  - z - Set of collocation/quadrature points
27  - w - Set of quadrature weights
28  - D - Derivative matrix
29  - h - Lagrange Interpolant
30  - I - Interpolation matrix
31  - g - Gauss
32  - k - Kronrod
33  - gr - Gauss-Radau
34  - gl - Gauss-Lobatto
35  - j - Jacobi
36  - m - point at minus 1 in Radau rules
37  - p - point at plus 1 in Radau rules
38 
39  -----------------------------------------------------------------------\n
40  MAIN ROUTINES\n
41  -----------------------------------------------------------------------\n
42 
43  Points and Weights:
44 
45  - zwgj Compute Gauss-Jacobi points and weights
46  - zwgrjm Compute Gauss-Radau-Jacobi points and weights (z=-1)
47  - zwgrjp Compute Gauss-Radau-Jacobi points and weights (z= 1)
48  - zwglj Compute Gauss-Lobatto-Jacobi points and weights
49  - zwgk Compute Gauss-Kronrod-Jacobi points and weights
50  - zwrk Compute Radau-Kronrod points and weights
51  - zwlk Compute Lobatto-Kronrod points and weights
52  - JacZeros Compute Gauss-Jacobi points and weights
53 
54  Derivative Matrices:
55 
56  - Dgj Compute Gauss-Jacobi derivative matrix
57  - Dgrjm Compute Gauss-Radau-Jacobi derivative matrix (z=-1)
58  - Dgrjp Compute Gauss-Radau-Jacobi derivative matrix (z= 1)
59  - Dglj Compute Gauss-Lobatto-Jacobi derivative matrix
60 
61  Lagrange Interpolants:
62 
63  - hgj Compute Gauss-Jacobi Lagrange interpolants
64  - hgrjm Compute Gauss-Radau-Jacobi Lagrange interpolants (z=-1)
65  - hgrjp Compute Gauss-Radau-Jacobi Lagrange interpolants (z= 1)
66  - hglj Compute Gauss-Lobatto-Jacobi Lagrange interpolants
67 
68  Interpolation Operators:
69 
70  - Imgj Compute interpolation operator gj->m
71  - Imgrjm Compute interpolation operator grj->m (z=-1)
72  - Imgrjp Compute interpolation operator grj->m (z= 1)
73  - Imglj Compute interpolation operator glj->m
74 
75  Polynomial Evaluation:
76 
77  - jacobfd Returns value and derivative of Jacobi poly. at point z
78  - jacobd Returns derivative of Jacobi poly. at point z (valid at z=-1,1)
79 
80  -----------------------------------------------------------------------\n
81  LOCAL ROUTINES\n
82  -----------------------------------------------------------------------\n
83 
84  - jacobz Returns Jacobi polynomial zeros
85  - gammaf Gamma function for integer values and halves
86  - RecCoeff Calculates the recurrence coefficients for orthogonal poly
87  - TriQL QL algorithm for symmetrix tridiagonal matrix
88  - JKMatrix Generates the Jacobi-Kronrod matrix
89 
90  ------------------------------------------------------------------------\n
91 
92  Useful references:
93 
94  - [1] Gabor Szego: Orthogonal Polynomials, American Mathematical Society,
95  Providence, Rhode Island, 1939.
96  - [2] Abramowitz \& Stegun: Handbook of Mathematical Functions,
97  Dover, New York, 1972.
98  - [3] Canuto, Hussaini, Quarteroni \& Zang: Spectral Methods in Fluid
99  Dynamics, Springer-Verlag, 1988.
100  - [4] Ghizzetti \& Ossicini: Quadrature Formulae, Academic Press, 1970.
101  - [5] Karniadakis \& Sherwin: Spectral/hp element methods for CFD, 1999
102 
103 
104  NOTES
105 
106  -# Legendre polynomial \f$ \alpha = \beta = 0 \f$
107  -# Chebychev polynomial \f$ \alpha = \beta = -0.5 \f$
108  -# All routines are double precision.
109  -# All array subscripts start from zero, i.e. vector[0..N-1]
110  */
111 
112  /*-----------------------------------------------------------------------
113  M A I N R O U T I N E S
114  -----------------------------------------------------------------------*/
115 
116  /* Points and weights */
117  LIB_UTILITIES_EXPORT void zwgj (double *, double *, const int , const double, const double);
118  LIB_UTILITIES_EXPORT void zwgrjm (double *, double *, const int , const double, const double);
119  LIB_UTILITIES_EXPORT void zwgrjp (double *, double *, const int , const double, const double);
120  LIB_UTILITIES_EXPORT void zwglj (double *, double *, const int , const double, const double);
121  LIB_UTILITIES_EXPORT void zwgk (double *, double *, const int , const double, const double);
122  LIB_UTILITIES_EXPORT void zwrk (double *, double *, const int , const double, const double);
123  LIB_UTILITIES_EXPORT void zwlk (double *, double *, const int , const double, const double);
124  LIB_UTILITIES_EXPORT void JacZeros(const int, double *, double*, const double,const double);
125 
126 
127  /* Derivative operators */
128  LIB_UTILITIES_EXPORT void Dgj (double *, const double *, const int, const double,
129  const double);
130  LIB_UTILITIES_EXPORT void Dgrjm (double *, const double *, const int, const double,
131  const double);
132  LIB_UTILITIES_EXPORT void Dgrjp (double *, const double *, const int, const double,
133  const double);
134  LIB_UTILITIES_EXPORT void Dglj (double *,const double *, const int, const double,
135  const double);
136 
137  /* Lagrangian interpolants */
138  LIB_UTILITIES_EXPORT double hgj (const int, const double, const double *, const int,
139  const double, const double);
140  LIB_UTILITIES_EXPORT double hgrjm (const int, const double, const double *, const int,
141  const double, const double);
142  LIB_UTILITIES_EXPORT double hgrjp (const int, const double, const double *, const int,
143  const double, const double);
144  LIB_UTILITIES_EXPORT double hglj (const int, const double, const double *, const int,
145  const double, const double);
146 
147  /* Interpolation operators */
148  LIB_UTILITIES_EXPORT void Imgj (double*, const double*, const double*, const int, const int,
149  const double, const double);
150  LIB_UTILITIES_EXPORT void Imgrjm(double*, const double*, const double*, const int, const int,
151  const double, const double);
152  LIB_UTILITIES_EXPORT void Imgrjp(double*, const double*, const double*, const int, const int,
153  const double, const double);
154  LIB_UTILITIES_EXPORT void Imglj (double*, const double*, const double*, const int, const int,
155  const double, const double);
156 
157  /* Polynomial functions */
158  LIB_UTILITIES_EXPORT void jacobfd (const int, const double *, double *, double *, const int ,
159  const double, const double);
160  LIB_UTILITIES_EXPORT void jacobd (const int, const double *, double *, const int ,
161  const double, const double);
162 
163  LIB_UTILITIES_EXPORT std::complex<Nektar::NekDouble> ImagBesselComp(const int,
164  std::complex<Nektar::NekDouble> );
165 
166 
167 } // end of namespace
168 
169 
170 #define H_PLYLIB
171 #endif /* END OF POLYLIB.H DECLARATIONS */
172 
173 
174 
175 
176 
177 
178 
179 
180 
void Dgj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated.
Definition: Polylib.cpp:937
void zwgrjp(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition: Polylib.cpp:239
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through.
Definition: Polylib.cpp:1345
void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points.
Definition: Polylib.cpp:1775
void zwglj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition: Polylib.cpp:311
double hgrjp(const int i, const double z, const double *zgrj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the.
Definition: Polylib.cpp:1505
void zwlk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Lobatto-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:711
void JacZeros(const int n, double *a, double *b, const double alpha, const double beta)
Zero and Weight determination through the eigenvalues and eigenvectors of a tridiagonal.
Definition: Polylib.cpp:2389
void jacobd(const int np, const double *z, double *polyd, const int n, const double alpha, const double beta)
Calculate the derivative of Jacobi polynomials.
Definition: Polylib.cpp:2151
void Imglj(double *im, const double *zglj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points.
Definition: Polylib.cpp:1837
#define LIB_UTILITIES_EXPORT
void Dgrjm(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated.
Definition: Polylib.cpp:1023
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:385
void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points.
Definition: Polylib.cpp:1715
void Dgrjp(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the.
Definition: Polylib.cpp:1121
void zwgrjm(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition: Polylib.cpp:163
The namespace associated with the the Polylib library (Polylib introduction)
Definition: Polylib.cpp:31
void zwrk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Radau-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:519
double hglj(const int i, const double z, const double *zglj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the.
Definition: Polylib.cpp:1589
std::complex< Nektar::NekDouble > ImagBesselComp(int n, std::complex< Nektar::NekDouble > y)
Calcualte the bessel function of the first kind with complex double input y. Taken from Numerical Rec...
Definition: Polylib.cpp:2999
double hgrjm(const int i, const double z, const double *zgrj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the.
Definition: Polylib.cpp:1421
void Dglj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the.
Definition: Polylib.cpp:1221
void Imgj(double *im, const double *zgj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Jacobi points to an.
Definition: Polylib.cpp:1657
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition: Polylib.cpp:1951
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:107