5 #include"Teuchos_LAPACK.hpp"
9 #ifndef __ORTHOGONAL_POLYNOMIALS__
10 #define __ORTHOGONAL_POLYNOMIALS__
29 void rec_jacobi( ROL::Ptr<Teuchos::LAPACK<int,Real> > lapack,
33 std::vector<Real> &b ) {
36 Real nu = (beta-alpha)/
double(alpha+beta+2.0);
37 Real mu = pow(2.0,alpha+beta+1.0)*tgamma(alpha+1.0)*tgamma(beta+1.0)/tgamma(alpha+beta+2.0);
39 Real sqdif = pow(beta,2)-pow(alpha,2);
46 for(
int n=1;n<N;++n) {
48 a[n] = sqdif/(nab*(nab+2));
51 b[1] = 4.0*(alpha+1.0)*(beta+1.0)/(pow(alpha+beta+2.0,2)*(alpha+beta+3.0));
54 for(
int n=2;n<N;++n) {
56 b[n] = 4.0*(n+alpha)*(n+beta)*n*(n+alpha+beta)/(nab*nab*(nab+1.0)*(nab-1.0));
73 const std::vector<Real> &b,
74 const std::vector<Real> &x,
75 std::vector<Real> &
V) {
76 const int n = a.size();
78 for(
int i=0;i<n;++i) {
84 for(
int j=1;j<n-1;++j) {
85 for(
int i=0;i<n;++i) {
86 V[i+(j+1)*n] = (x[i] - a[j])*
V[i+j*n] - b[j]*
V[i+(j-1)*n];
104 void gauss( ROL::Ptr<Teuchos::LAPACK<int,Real> > lapack,
105 const std::vector<Real> &a,
106 const std::vector<Real> &b,
107 std::vector<Real> &x,
108 std::vector<Real> &w ) {
110 const int N = a.size();
112 const char COMPZ =
'I';
114 std::vector<Real> D(N,0.0);
115 std::vector<Real> E(N,0.0);
116 std::vector<Real> WORK(4*N,0.0);
119 std::vector<Real> Z(N*N,0);
123 for(
int i=0;i<N-1;++i) {
127 lapack->STEQR(COMPZ,N,&D[0],&E[0],&Z[0],LDZ,&WORK[0],&INFO);
129 for(
int i=0;i<N;++i){
131 w[i] = b[0]*pow(Z[N*i],2);
151 void rec_lobatto( ROL::Ptr<Teuchos::LAPACK<int,Real> >
const lapack,
154 std::vector<Real> &a,
155 std::vector<Real> &b ) {
156 const int N = a.size()-1;
158 std::vector<Real> amod(N,0.0);
159 std::vector<Real> bmod(N-1,0.0);
160 std::vector<Real> en(N,0.0);
161 std::vector<Real> g(N,0.0);
166 for(
int i=0;i<N-1;++i) {
167 bmod[i] = sqrt(b[i+1]);
170 for(
int i=0;i<N;++i) {
174 trisolve(lapack,bmod,amod,bmod,en,g);
177 for(
int i=0;i<N;++i) {
181 trisolve(lapack,bmod,amod,bmod,en,g);
184 a[N] = (g1*xl2-g2*xl1)/(g1-g2);
185 b[N] = (xl2-xl1)/(g1-g2);