33 #ifndef OOMPH_ORTHPOLY_HEADER 34 #define OOMPH_ORTHPOLY_HEADER 38 #include <oomph-lib-config.h> 57 const double eps = 1.0e-15;
62 inline double legendre(
const unsigned &p,
const double &x)
67 else if(p==1)
return x;
72 double L0 = 1.0, L1 = 1.0, L2 = x;
74 for(
unsigned n=1;n<p;n++)
79 L2 = ((2*n + 1)*x*L1 - n*L0)/(n+1);
100 double L0 = 1.0, L1 = 1.0, L2 = x;
102 for(
unsigned n=1;n<p;n++){
107 L2 = ((2*n+1)*x*L1 - n*L0)/(n+1);
118 inline double dlegendre(
const unsigned &p,
const double &x)
120 double dL1 = 1.0, dL2 = 3*x, dL3=0.0;
122 else if(p==1)
return dL1;
123 else if(p==2)
return dL2;
125 for(
unsigned n=2;n<p;n++){
126 dL3 = 1.0/n*((2.0*n+1.0)*x*dL2-(n+1.0)*dL1);
139 double ddL2 = 3.0, ddL3 = 15*x, ddL4=0.0;
141 else if(p==1)
return 0.0;
142 else if(p==2)
return ddL2;
143 else if(p==3)
return ddL3;
145 for(
unsigned i=3;
i<p;
i++){
146 ddL4 = 1.0/(
i-1.0)*((2.0*
i+1.0)*x*ddL3-(
i+2.0)*ddL2);
155 inline double jacobi(
const int &alpha,
const int &beta,
const unsigned &p,
159 double P1 = 0.5*(alpha-beta+(alpha+beta+2.0)*x);
162 else if(p==1)
return P1;
164 for(
unsigned n=1;n<p;n++){
165 double a1n = 2*(n+1)*(n+alpha+beta+1)*(2*n+alpha+beta);
166 double a2n = (2*n+alpha+beta+1)*(alpha*alpha-beta*beta);
167 double a3n = (2*n+alpha+beta)*(2*n+alpha+beta+1)*(2*n+alpha+beta+2);
168 double a4n = 2*(n+alpha)*(n+beta)*(2*n+alpha+beta+2);
169 P2 = ( (a2n+a3n*x)*P1 - a4n*P0 ) / a1n;
178 inline void jacobi(
const int &alpha,
const int &beta,
const unsigned &p,
187 polys[1] = 0.5*(alpha-beta+(alpha+beta+2.0)*x);
191 double P2 = 0.5*(alpha-beta+(alpha+beta+2.0)*x);
193 for(
unsigned n=1;n<p;n++)
199 double a1n = 2*(n+1)*(n+alpha+beta+1)*(2*n+alpha+beta);
200 double a2n = (2*n+alpha+beta+1)*(alpha*alpha-beta*beta);
201 double a3n = (2*n+alpha+beta)*(2*n+alpha+beta+1)*(2*n+alpha+beta+2);
202 double a4n = 2*(n+alpha)*(n+beta)*(2*n+alpha+beta+2);
204 P2 = ( (a2n+a3n*x)*P1 - a4n*P0 ) / a1n;
double dlegendre(const unsigned &p, const double &x)
Calculates first derivative of Legendre polynomial of degree p at x using three term recursive formul...
double legendre(const unsigned &p, const double &x)
Calculates Legendre polynomial of degree p at x using the three term recurrence relation ...
void gll_nodes(const unsigned &Nnode, Vector< double > &x)
Calculates the Gauss Lobatto Legendre abscissas for degree p = NNode-1.
double ddlegendre(const unsigned &p, const double &x)
Calculates second derivative of Legendre polynomial of degree p at x using three term recursive formu...
void gl_nodes(const unsigned &Nnode, Vector< double > &x)
void legendre_vector(const unsigned &p, const double &x, Vector< double > &polys)
Calculates Legendre polynomial of degree p at x using three term recursive formula. Returns all polynomials up to order p in the vector.
double jacobi(const int &alpha, const int &beta, const unsigned &p, const double &x)
Calculate the Jacobi polnomials.