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