00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040 #include <iostream>
00041
00042 using namespace std ;
00043
00044 #include <stdlib.h>
00045 #include <math.h>
00046 #include <assert.h>
00047
00048 #include "ortho_poly.h"
00049 #include "grid.h"
00050
00051 void legendre_nodes_weight_GL(int n, double* coloc, double* weight,
00052 double prec = 1.e-12, int itemax = 100) ;
00053
00054
00055
00056
00057
00058
00059 Legendre_poly::Legendre_poly(int ni) : Ortho_poly(ni) {}
00060
00061
00062
00063
00064
00065 Legendre_poly::Legendre_poly(const Legendre_poly& pp) : Ortho_poly(pp) {}
00066
00067
00068
00069
00070
00071
00072 Legendre_poly::~Legendre_poly() {}
00073
00074
00075
00076
00077
00078
00079 void Legendre_poly::operator=(const Legendre_poly& ) {
00080
00081 cerr << "Legendre_poly::operator= not implemented !" << endl ;
00082 abort() ;
00083
00084 }
00085
00086
00087
00088
00089
00090 ostream& operator<<(ostream& ost, const Legendre_poly& pp) {
00091
00092 ost << "Basis of Legendre polynomials up to degree " << pp.n() << endl ;
00093
00094 return ost ;
00095 }
00096
00097
00098
00099
00100
00101 double Legendre_poly::weight(double ) const {
00102
00103 return 1. ;
00104
00105 }
00106
00107
00108
00109
00110
00111
00112 double Legendre_poly::operator()(int i, double x) const {
00113
00114 assert( i >= 0 ) ;
00115 assert( i <= nn ) ;
00116
00117 if (i==0) return 1. ;
00118
00119 if (i==1) return x ;
00120
00121 double tjm2 = 1. ;
00122 double tjm1 = x ;
00123
00124 for (int j=2; j<=i; j++) {
00125 double tj = ( (2*j-1)*x* tjm1 - (j-1)*tjm2 ) / double(j) ;
00126 tjm2 = tjm1 ;
00127 tjm1 = tj ;
00128 }
00129
00130 return tjm1 ;
00131
00132 }
00133
00134
00135
00136
00137
00138
00139 const Grid& Legendre_poly::gauss_nodes() const {
00140
00141 if (p_gauss_nodes == 0x0) {
00142
00143 p_gauss_nodes = new Grid_Legendre_Gauss(nn+1) ;
00144 }
00145
00146 return *p_gauss_nodes ;
00147 }
00148
00149
00150 double Legendre_poly::gauss_weight(int i) const {
00151
00152 if (p_gauss_weights == 0x0) {
00153
00154 cerr << "Legendre_poly::gauss_weight : not implemented yet !"
00155 << endl ;
00156 abort() ;
00157 }
00158
00159 return p_gauss_weights[i] ;
00160
00161 }
00162
00163
00164 double Legendre_poly::gauss_gamma(int i) const {
00165
00166 if (p_gauss_gamma == 0x0) {
00167
00168 p_gauss_gamma = new double[nn+1] ;
00169
00170 for (int j=0; j<=nn; j++)
00171 p_gauss_gamma[j] = 1. / (double(j) + 0.5) ;
00172 }
00173
00174 return p_gauss_gamma[i] ;
00175
00176 }
00177
00178
00179 const Grid& Legendre_poly::gauss_lobatto_nodes() const {
00180
00181 if (p_gauss_lobatto_nodes == 0x0) {
00182
00183 p_gauss_lobatto_nodes = new Grid_Legendre_GL(nn+1) ;
00184 }
00185
00186 return *p_gauss_lobatto_nodes ;
00187 }
00188
00189
00190 double Legendre_poly::gauss_lobatto_weight(int i) const {
00191
00192 if (p_gauss_lobatto_weights == 0x0) {
00193
00194 p_gauss_lobatto_weights = new double[nn+1] ;
00195
00196 double* xx_work = new double[nn+1] ;
00197
00198 double precis = 1.e-12 ;
00199 int nitermax = 200 ;
00200
00201 legendre_nodes_weight_GL(nn+1, xx_work, p_gauss_lobatto_weights,
00202 precis, nitermax) ;
00203
00204 delete [] xx_work ;
00205 }
00206
00207 return p_gauss_lobatto_weights[i] ;
00208
00209 }
00210
00211
00212 double Legendre_poly::gauss_lobatto_gamma(int i) const {
00213
00214 if (p_gauss_lobatto_gamma == 0x0) {
00215
00216 p_gauss_lobatto_gamma = new double[nn+1] ;
00217
00218 for (int j=0; j<nn; j++)
00219 p_gauss_lobatto_gamma[j] = 1. / (double(j) + 0.5) ;
00220
00221 p_gauss_lobatto_gamma[nn] = 2. / double(nn) ;
00222
00223 }
00224
00225 return p_gauss_lobatto_gamma[i] ;
00226
00227 }
00228
00229
00230
00231
00232
00233
00234 void Legendre_poly::coef_projection(double (*f)(double), double* cf) const {
00235
00236
00237
00238
00239
00240 int n_large = 128 ;
00241
00242 if (nn > n_large) {
00243 cerr << "Legendre_poly::coef_projection : nn > n_large !"
00244 << endl << " nn = " << nn << " n_large = " << n_large << endl ;
00245 abort() ;
00246 }
00247
00248 Legendre_poly leg_large(n_large) ;
00249
00250 double* cf_large = new double[n_large + 1] ;
00251
00252 leg_large.coef_interpolant_GL(f, cf_large) ;
00253
00254 for (int i=0; i<=nn; i++) cf[i] = cf_large[i] ;
00255
00256 delete [] cf_large ;
00257
00258 }
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270