00001 /* 00002 * Copyright (c) 1999-2001 Eric Gourgoulhon 00003 * 00004 * This file is part of LORENE. 00005 * 00006 * LORENE is free software; you can redistribute it and/or modify 00007 * it under the terms of the GNU General Public License as published by 00008 * the Free Software Foundation; either version 2 of the License, or 00009 * (at your option) any later version. 00010 * 00011 * LORENE is distributed in the hope that it will be useful, 00012 * but WITHOUT ANY WARRANTY; without even the implied warranty of 00013 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00014 * GNU General Public License for more details. 00015 * 00016 * You should have received a copy of the GNU General Public License 00017 * along with LORENE; if not, write to the Free Software 00018 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00019 * 00020 */ 00021 00022 00023 char circhebpimp_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/circhebpimp.C,v 1.2 2004/12/27 14:27:28 j_novak Exp $" ; 00024 00025 00026 /* 00027 * Transformation de Tchebyshev inverse T_{2k}/T_{2k+1} (suivant la parite de 00028 * l'indice m en phi) sur le troisieme indice 00029 * (indice correspondant a r) d'un tableau 3-D decrivant une fonction symetrique 00030 * par rapport au plan equatorial z = 0 et sans aucune autre symetrie, 00031 * cad que l'on a effectue 00032 * 1/ un developpement en polynomes de Tchebyshev pairs pour m pair 00033 * 2/ un developpement en polynomes de Tchebyshev impairs pour m impair 00034 * 00035 * 00036 * Utilise la bibliotheque fftw. 00037 * 00038 * Entree: 00039 * ------- 00040 * int* deg : tableau du nombre effectif de degres de liberte dans chacune 00041 * des 3 dimensions: le nombre de points de collocation 00042 * en r est nr = deg[2] et doit etre de la forme 00043 * nr = 2*p + 1 00044 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois 00045 * dimensions. 00046 * On doit avoir dimc[2] >= deg[2] = nr. 00047 * 00048 * double* cf : tableau des coefficients c_i de la fonction definis 00049 * comme suit (a theta et phi fixes) 00050 * 00051 * -- pour m pair (i.e. j = 0, 1, 4, 5, 8, 9, ...) : 00052 * 00053 * f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) , 00054 * 00055 * ou T_{2i}(x) designe le polynome de Tchebyshev de 00056 * degre 2i. 00057 * 00058 * -- pour m impair (i.e. j = 2, 3, 6, 7, 10, 11, ...) : 00059 * 00060 * f(x) = som_{i=0}^{nr-2} c_i T_{2i+1}(x) , 00061 * 00062 * ou T_{2i+1}(x) designe le polynome de Tchebyshev de 00063 * degre 2i+1. 00064 * 00065 * Les coefficients c_i (0 <= i <= nr-1) doivent etre stokes 00066 * dans le tableau cf comme suit 00067 * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ] 00068 * ou j et k sont les indices correspondant a phi et theta 00069 * respectivement. 00070 * L'espace memoire correspondant a ce pointeur doit etre 00071 * dimc[0]*dimc[1]*dimc[2] et doit etre alloue avant l'appel a 00072 * la routine. 00073 * 00074 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois 00075 * dimensions. 00076 * On doit avoir dimf[2] >= deg[2] = nr. 00077 * 00078 * Sortie: 00079 * ------- 00080 * double* ff : tableau des valeurs de la fonction aux nr points de 00081 * de collocation 00082 * 00083 * x_i = sin( pi/2 i/(nr-1) ) 0 <= i <= nr-1 00084 * 00085 * Les valeurs de la fonction sont stokees dans le 00086 * tableau ff comme suit 00087 * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ] 00088 * ou j et k sont les indices correspondant a phi et theta 00089 * respectivement. 00090 * L'espace memoire correspondant a ce pointeur doit etre 00091 * dimf[0]*dimf[1]*dimf[2] et doit avoir ete alloue avant 00092 * l'appel a la routine. 00093 * 00094 * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un 00095 * seul tableau, qui constitue une entree/sortie. 00096 */ 00097 00098 /* 00099 * $Id: circhebpimp.C,v 1.2 2004/12/27 14:27:28 j_novak Exp $ 00100 * $Log: circhebpimp.C,v $ 00101 * Revision 1.2 2004/12/27 14:27:28 j_novak 00102 * Added forgotten "delete [] t1" 00103 * 00104 * Revision 1.1 2004/12/21 17:06:03 j_novak 00105 * Added all files for using fftw3. 00106 * 00107 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon 00108 * Suppressed the directive #include <malloc.h> for malloc is defined 00109 * in <stdlib.h> 00110 * 00111 * Revision 1.3 2002/10/16 14:36:53 j_novak 00112 * Reorganization of #include instructions of standard C++, in order to 00113 * use experimental version 3 of gcc. 00114 * 00115 * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon 00116 * Modification of declaration of Fortran 77 prototypes for 00117 * a better portability (in particular on IBM AIX systems): 00118 * All Fortran subroutine names are now written F77_* and are 00119 * defined in the new file C++/Include/proto_f77.h. 00120 * 00121 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon 00122 * LORENE 00123 * 00124 * Revision 2.0 1999/02/22 15:43:10 hyc 00125 * *** empty log message *** 00126 * 00127 * 00128 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/circhebpimp.C,v 1.2 2004/12/27 14:27:28 j_novak Exp $ 00129 * 00130 */ 00131 00132 00133 // headers du C 00134 #include <assert.h> 00135 #include <stdlib.h> 00136 #include <fftw3.h> 00137 00138 //Lorene prototypes 00139 #include "tbl.h" 00140 00141 // Prototypage des sous-routines utilisees: 00142 fftw_plan back_fft(int, Tbl*&) ; 00143 double* cheb_ini(const int) ; 00144 double* chebimp_ini(const int ) ; 00145 //***************************************************************************** 00146 00147 void circhebpimp(const int* deg, const int* dimc, double* cf, 00148 const int* dimf, double* ff) 00149 00150 { 00151 int i, j, k ; 00152 00153 // Dimensions des tableaux ff et cf : 00154 int n1f = dimf[0] ; 00155 int n2f = dimf[1] ; 00156 int n3f = dimf[2] ; 00157 int n1c = dimc[0] ; 00158 int n2c = dimc[1] ; 00159 int n3c = dimc[2] ; 00160 00161 // Nombres de degres de liberte en r : 00162 int nr = deg[2] ; 00163 00164 // Tests de dimension: 00165 if (nr > n3c) { 00166 cout << "circhebpimp: nr > n3c : nr = " << nr << " , n3c = " 00167 << n3c << endl ; 00168 abort () ; 00169 exit(-1) ; 00170 } 00171 if (nr > n3f) { 00172 cout << "circhebpimp: nr > n3f : nr = " << nr << " , n3f = " 00173 << n3f << endl ; 00174 abort () ; 00175 exit(-1) ; 00176 } 00177 if (n1c > n1f) { 00178 cout << "circhebpimp: n1c > n1f : n1c = " << n1c << " , n1f = " 00179 << n1f << endl ; 00180 abort () ; 00181 exit(-1) ; 00182 } 00183 if (n2c > n2f) { 00184 cout << "circhebpimp: n2c > n2f : n2c = " << n2c << " , n2f = " 00185 << n2f << endl ; 00186 abort () ; 00187 exit(-1) ; 00188 } 00189 00190 // Nombre de points pour la FFT: 00191 int nm1 = nr - 1; 00192 int nm1s2 = nm1 / 2; 00193 00194 // Recherche des tables pour la FFT: 00195 Tbl* pg = 0x0 ; 00196 fftw_plan p = back_fft(nm1, pg) ; 00197 Tbl& g = *pg ; 00198 double* t1 = new double[nr] ; 00199 00200 // Recherche de la table des sin(psi) : 00201 double* sinp = cheb_ini(nr); 00202 00203 // Recherche de la table des points de collocations x_k : 00204 double* x = chebimp_ini(nr); 00205 00206 // boucle sur phi et theta 00207 00208 int n2n3f = n2f * n3f ; 00209 int n2n3c = n2c * n3c ; 00210 00211 //======================================================================= 00212 // Cas m pair 00213 //======================================================================= 00214 00215 j = 0 ; 00216 00217 while (j<n1f-1) { //le dernier coef en phi n'est pas considere 00218 // (car nul) 00219 00220 //-------------------------------------------------------------------- 00221 // partie cos(m phi) avec m pair : developpement en T_{2i}(x) 00222 //-------------------------------------------------------------------- 00223 00224 for (k=0; k<n2c; k++) { 00225 00226 int i0 = n2n3c * j + n3c * k ; // indice de depart 00227 double* cf0 = cf + i0 ; // tableau des donnees a transformer 00228 00229 i0 = n2n3f * j + n3f * k ; // indice de depart 00230 double* ff0 = ff + i0 ; // tableau resultat 00231 00232 /* 00233 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi] 00234 * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)). 00235 */ 00236 00237 // Calcul des coefficients de Fourier de la fonction 00238 // G(psi) = F+(psi) + F_(psi) sin(psi) 00239 // en fonction des coefficients de Tchebyshev de f: 00240 00241 // Coefficients impairs de G 00242 //-------------------------- 00243 00244 double c1 = cf0[1] ; 00245 00246 double som = 0; 00247 ff0[1] = 0 ; 00248 for ( i = 3; i < nr; i += 2 ) { 00249 ff0[i] = cf0[i] - c1 ; 00250 som += ff0[i] ; 00251 } 00252 00253 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00254 double fmoins0 = nm1s2 * c1 + som ; 00255 00256 // Coef. impairs de G 00257 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw 00258 // donnait exactement les coef. des sinus, ce facteur serait -0.5. 00259 for ( i = 3; i < nr; i += 2 ) { 00260 g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ; 00261 } 00262 00263 00264 // Coefficients pairs de G 00265 //------------------------ 00266 // Ces coefficients sont egaux aux coefficients pairs du developpement de 00267 // f. 00268 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw 00269 // donnait exactement les coef. des cosinus, ce facteur serait 1. 00270 00271 g.set(0) = cf0[0] ; 00272 for (i=1; i<nm1s2; i++) g.set(i) = 0.5 * cf0[2*i] ; 00273 g.set(nm1s2) = cf0[nm1] ; 00274 00275 // Transformation de Fourier inverse de G 00276 //--------------------------------------- 00277 00278 // FFT inverse 00279 fftw_execute(p) ; 00280 00281 // Valeurs de f deduites de celles de G 00282 //------------------------------------- 00283 00284 for ( i = 1; i < nm1s2 ; i++ ) { 00285 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2: 00286 int isym = nm1 - i ; 00287 // ... indice (dans le tableau ff0) du point x correspondant a psi 00288 int ix = nm1 - i ; 00289 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi) 00290 int ixsym = nm1 - isym ; 00291 00292 double fp = .5 * ( g(i) + g(isym) ) ; 00293 double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ; 00294 00295 ff0[ix] = fp + fm ; 00296 ff0[ixsym] = fp - fm ; 00297 } 00298 00299 //... cas particuliers: 00300 ff0[0] = g(0) - fmoins0 ; 00301 ff0[nm1] = g(0) + fmoins0 ; 00302 ff0[nm1s2] = g(nm1s2) ; 00303 00304 } // fin de la boucle sur theta 00305 00306 //-------------------------------------------------------------------- 00307 // partie sin(m phi) avec m pair : developpement en T_{2i}(x) 00308 //-------------------------------------------------------------------- 00309 00310 j++ ; 00311 00312 if ( (j != 1) && (j != n1f-1) ) { 00313 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont 00314 // pas nuls 00315 00316 for (k=0; k<n2c; k++) { 00317 00318 int i0 = n2n3c * j + n3c * k ; // indice de depart 00319 double* cf0 = cf + i0 ; // tableau des donnees a transformer 00320 00321 i0 = n2n3f * j + n3f * k ; // indice de depart 00322 double* ff0 = ff + i0 ; // tableau resultat 00323 00324 /* 00325 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi] 00326 * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)). 00327 */ 00328 00329 // Calcul des coefficients de Fourier de la fonction 00330 // G(psi) = F+(psi) + F_(psi) sin(psi) 00331 // en fonction des coefficients de Tchebyshev de f: 00332 00333 // Coefficients impairs de G 00334 //-------------------------- 00335 00336 double c1 = cf0[1] ; 00337 00338 double som = 0; 00339 ff0[1] = 0 ; 00340 for ( i = 3; i < nr; i += 2 ) { 00341 ff0[i] = cf0[i] - c1 ; 00342 som += ff0[i] ; 00343 } 00344 00345 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00346 double fmoins0 = nm1s2 * c1 + som ; 00347 00348 // Coef. impairs de G 00349 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw 00350 // donnait exactement les coef. des sinus, ce facteur serait -0.5. 00351 for ( i = 3; i < nr; i += 2 ) { 00352 g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ; 00353 } 00354 00355 00356 // Coefficients pairs de G 00357 //------------------------ 00358 // Ces coefficients sont egaux aux coefficients pairs du developpement de 00359 // f. 00360 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw 00361 // donnait exactement les coef. des cosinus, ce facteur serait 1. 00362 00363 g.set(0) = cf0[0] ; 00364 for (i=1; i<nm1s2; i++) g.set(i) = 0.5 * cf0[2*i] ; 00365 g.set(nm1s2) = cf0[nm1] ; 00366 00367 // Transformation de Fourier inverse de G 00368 //--------------------------------------- 00369 00370 // FFT inverse 00371 fftw_execute(p) ; 00372 00373 // Valeurs de f deduites de celles de G 00374 //------------------------------------- 00375 00376 for ( i = 1; i < nm1s2 ; i++ ) { 00377 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2: 00378 int isym = nm1 - i ; 00379 // ... indice (dans le tableau ff0) du point x correspondant a psi 00380 int ix = nm1 - i ; 00381 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi) 00382 int ixsym = nm1 - isym ; 00383 00384 double fp = .5 * ( g(i) + g(isym) ) ; 00385 double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ; 00386 00387 ff0[ix] = fp + fm ; 00388 ff0[ixsym] = fp - fm ; 00389 } 00390 00391 //... cas particuliers: 00392 ff0[0] = g(0) - fmoins0 ; 00393 ff0[nm1] = g(0) + fmoins0 ; 00394 ff0[nm1s2] = g(nm1s2) ; 00395 00396 } // fin de la boucle sur theta 00397 00398 } // fin du cas ou le calcul etait necessaire (i.e. ou les 00399 // coef en phi n'etaient pas nuls) 00400 00401 // On passe au cas m pair suivant: 00402 j+=3 ; 00403 00404 } // fin de la boucle sur les cas m pair 00405 00406 if (n1f<=3) { // cas m=0 seulement (symetrie axiale) 00407 delete [] t1 ; 00408 return ; 00409 } 00410 00411 //======================================================================= 00412 // Cas m impair 00413 //======================================================================= 00414 00415 j = 2 ; 00416 00417 while (j<n1f-1) { //le dernier coef en phi n'est pas considere 00418 // (car nul) 00419 00420 //------------------------------------------------------------------------ 00421 // partie cos(m phi) avec m impair : developpement en T_{2i+1}(x) 00422 //------------------------------------------------------------------------ 00423 00424 for (k=0; k<n2c; k++) { 00425 00426 int i0 = n2n3c * j + n3c * k ; // indice de depart 00427 double* cf0 = cf + i0 ; // tableau des donnees a transformer 00428 00429 i0 = n2n3f * j + n3f * k ; // indice de depart 00430 double* ff0 = ff + i0 ; // tableau resultat 00431 00432 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction 00433 // h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le 00434 // tableau t1 : 00435 00436 t1[0] = .5 * cf0[0] ; 00437 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ; 00438 t1[nm1] = .5 * cf0[nr-2] ; 00439 00440 /* 00441 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi] 00442 * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)). 00443 */ 00444 00445 // Calcul des coefficients de Fourier de la fonction 00446 // G(psi) = F+(psi) + F_(psi) sin(psi) 00447 // en fonction des coefficients de Tchebyshev de f: 00448 00449 // Coefficients impairs de G 00450 //-------------------------- 00451 00452 double c1 = t1[1] ; 00453 00454 double som = 0; 00455 ff0[1] = 0 ; 00456 for ( i = 3; i < nr; i += 2 ) { 00457 ff0[i] = t1[i] - c1 ; 00458 som += ff0[i] ; 00459 } 00460 00461 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00462 double fmoins0 = nm1s2 * c1 + som ; 00463 00464 // Coef. impairs de G 00465 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw 00466 // donnait exactement les coef. des sinus, ce facteur serait -0.5. 00467 for ( i = 3; i < nr; i += 2 ) { 00468 g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ; 00469 } 00470 00471 00472 // Coefficients pairs de G 00473 //------------------------ 00474 // Ces coefficients sont egaux aux coefficients pairs du developpement de 00475 // f. 00476 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw 00477 // donnait exactement les coef. des cosinus, ce facteur serait 1. 00478 00479 g.set(0) = t1[0] ; 00480 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ; 00481 g.set(nm1s2) = t1[nm1] ; 00482 00483 // Transformation de Fourier inverse de G 00484 //--------------------------------------- 00485 00486 // FFT inverse 00487 fftw_execute(p) ; 00488 00489 // Valeurs de f deduites de celles de G 00490 //------------------------------------- 00491 00492 for ( i = 1; i < nm1s2 ; i++ ) { 00493 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2: 00494 int isym = nm1 - i ; 00495 // ... indice (dans le tableau ff0) du point x correspondant a psi 00496 int ix = nm1 - i ; 00497 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi) 00498 int ixsym = nm1 - isym ; 00499 00500 double fp = .5 * ( g(i) + g(isym) ) ; 00501 double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ; 00502 00503 ff0[ix] = ( fp + fm ) / x[ix]; 00504 ff0[ixsym] = ( fp - fm ) / x[ixsym] ; 00505 } 00506 00507 //... cas particuliers: 00508 ff0[0] = 0 ; 00509 ff0[nm1] = g(0) + fmoins0 ; 00510 ff0[nm1s2] = g(nm1s2) / x[nm1s2] ; 00511 } // fin de la boucle sur theta 00512 00513 //------------------------------------------------------------------------ 00514 // partie sin(m phi) avec m impair : developpement en T_{2i+1}(x) 00515 //------------------------------------------------------------------------ 00516 00517 j++ ; 00518 00519 if ( j != n1f-1 ) { 00520 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont 00521 // pas nuls 00522 00523 for (k=0; k<n2c; k++) { 00524 00525 int i0 = n2n3c * j + n3c * k ; // indice de depart 00526 double* cf0 = cf + i0 ; // tableau des donnees a transformer 00527 00528 i0 = n2n3f * j + n3f * k ; // indice de depart 00529 double* ff0 = ff + i0 ; // tableau resultat 00530 00531 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction 00532 // h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le 00533 // tableau t1 : 00534 t1[0] = .5 * cf0[0] ; 00535 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ; 00536 t1[nm1] = .5 * cf0[nr-2] ; 00537 00538 /* 00539 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi] 00540 * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)). 00541 */ 00542 00543 // Calcul des coefficients de Fourier de la fonction 00544 // G(psi) = F+(psi) + F_(psi) sin(psi) 00545 // en fonction des coefficients de Tchebyshev de f: 00546 00547 // Coefficients impairs de G 00548 //-------------------------- 00549 00550 double c1 = t1[1] ; 00551 00552 double som = 0; 00553 ff0[1] = 0 ; 00554 for ( i = 3; i < nr; i += 2 ) { 00555 ff0[i] = t1[i] - c1 ; 00556 som += ff0[i] ; 00557 } 00558 00559 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00560 double fmoins0 = nm1s2 * c1 + som ; 00561 00562 // Coef. impairs de G 00563 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw 00564 // donnait exactement les coef. des sinus, ce facteur serait -0.5. 00565 for ( i = 3; i < nr; i += 2 ) { 00566 g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ; 00567 } 00568 00569 00570 // Coefficients pairs de G 00571 //------------------------ 00572 // Ces coefficients sont egaux aux coefficients pairs du developpement de 00573 // f. 00574 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw 00575 // donnait exactement les coef. des cosinus, ce facteur serait 1. 00576 00577 g.set(0) = t1[0] ; 00578 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ; 00579 g.set(nm1s2) = t1[nm1] ; 00580 00581 // Transformation de Fourier inverse de G 00582 //--------------------------------------- 00583 00584 // FFT inverse 00585 fftw_execute(p) ; 00586 00587 // Valeurs de f deduites de celles de G 00588 //------------------------------------- 00589 00590 for ( i = 1; i < nm1s2 ; i++ ) { 00591 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2: 00592 int isym = nm1 - i ; 00593 // ... indice (dans le tableau ff0) du point x correspondant a psi 00594 int ix = nm1 - i ; 00595 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi) 00596 int ixsym = nm1 - isym ; 00597 00598 double fp = .5 * ( g(i) + g(isym) ) ; 00599 double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ; 00600 00601 ff0[ix] = ( fp + fm ) / x[ix]; 00602 ff0[ixsym] = ( fp - fm ) / x[ixsym] ; 00603 } 00604 00605 //... cas particuliers: 00606 ff0[0] = 0 ; 00607 ff0[nm1] = g(0) + fmoins0 ; 00608 ff0[nm1s2] = g(nm1s2) / x[nm1s2] ; 00609 00610 } // fin de la boucle sur theta 00611 00612 } // fin du cas ou le calcul etait necessaire (i.e. ou les 00613 // coef en phi n'etaient pas nuls) 00614 00615 // On passe au cas m impair suivant: 00616 j+=3 ; 00617 00618 } // fin de la boucle sur les cas m impair 00619 00620 delete [] t1 ; 00621 }
1.4.6