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