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 cftcossins_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossins.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $" ; 00024 00025 /* 00026 * Transformation en sin(l*theta) ou cos(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 sans symetrie par rapport 00029 * au plan z=0. 00030 * Utilise la bibliotheque 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 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-1} c_l sin(l theta ) . 00069 * pour m impair: 00070 * f(theta) = som_{l=0}^{nt-1} c_l cos( 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_0 = c_{nt-1} = 0. 00081 * Pour m impair, c_{nt-1} = 0. 00082 * 00083 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un 00084 * seul tableau, qui constitue une entree/sortie. 00085 * 00086 */ 00087 00088 /* 00089 * $Id: cftcossins.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00090 * $Log: cftcossins.C,v $ 00091 * Revision 1.1 2004/12/21 17:06:02 j_novak 00092 * Added all files for using fftw3. 00093 * 00094 * Revision 1.1 2004/11/23 15:13:50 m_forot 00095 * Added the bases for the cases without any equatorial symmetry 00096 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I). 00097 * 00098 * 00099 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossins.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00100 * 00101 */ 00102 00103 // headers du C 00104 #include <stdlib.h> 00105 #include <fftw3.h> 00106 00108 #include "tbl.h" 00109 00110 // Prototypage des sous-routines utilisees: 00111 fftw_plan prepare_fft(int, Tbl*&) ; 00112 double* cheb_ini(const int) ; 00113 //***************************************************************************** 00114 00115 void cftcossins(const int* deg, const int* dimf, double* ff, const int* dimc, 00116 double* cf) 00117 { 00118 00119 int i, j, k ; 00120 00121 // Dimensions des tableaux ff et cf : 00122 int n1f = dimf[0] ; 00123 int n2f = dimf[1] ; 00124 int n3f = dimf[2] ; 00125 int n1c = dimc[0] ; 00126 int n2c = dimc[1] ; 00127 int n3c = dimc[2] ; 00128 00129 // Nombre de degres de liberte en theta : 00130 int nt = deg[1] ; 00131 00132 // Tests de dimension: 00133 if (nt > n2f) { 00134 cout << "cftcossins: nt > n2f : nt = " << nt << " , n2f = " 00135 << n2f << endl ; 00136 abort () ; 00137 exit(-1) ; 00138 } 00139 if (nt > n2c) { 00140 cout << "cftcossins: nt > n2c : nt = " << nt << " , n2c = " 00141 << n2c << endl ; 00142 abort () ; 00143 exit(-1) ; 00144 } 00145 if (n1f > n1c) { 00146 cout << "cftcossins: n1f > n1c : n1f = " << n1f << " , n1c = " 00147 << n1c << endl ; 00148 abort () ; 00149 exit(-1) ; 00150 } 00151 if (n3f > n3c) { 00152 cout << "cftcossins: n3f > n3c : n3f = " << n3f << " , n3c = " 00153 << n3c << endl ; 00154 abort () ; 00155 exit(-1) ; 00156 } 00157 00158 // Nombre de points pour la FFT: 00159 int nm1 = nt - 1; 00160 int nm1s2 = nm1 / 2; 00161 00162 // Recherche des tables pour la FFT: 00163 Tbl* pg = 0x0 ; 00164 fftw_plan p = prepare_fft(nm1, pg) ; 00165 Tbl& g = *pg ; 00166 00167 // Recherche de la table des sin(psi) : 00168 double* sinp = cheb_ini(nt); 00169 00170 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1 00171 // et 0 a dimf[2]) 00172 00173 int n2n3f = n2f * n3f ; 00174 int n2n3c = n2c * n3c ; 00175 00176 //======================================================================= 00177 // Cas m pair 00178 //======================================================================= 00179 00180 j = 0 ; 00181 00182 while (j<n1f-1) { //le dernier coef en phi n'est pas considere 00183 // (car nul) 00184 00185 //-------------------------------------------------------------------- 00186 // partie cos(m phi) avec m pair : transformation en sin(l) theta) 00187 //-------------------------------------------------------------------- 00188 00189 for (k=0; k<n3f; k++) { 00190 00191 int i0 = n2n3f * j + k ; // indice de depart 00192 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00193 00194 i0 = n2n3c * j + k ; // indice de depart 00195 double* cf0 = cf + i0 ; // tableau resultat 00196 00197 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta) 00198 //--------------------------------------------- 00199 00200 for ( i = 1; i < nm1s2 ; i++ ) { 00201 int isym = nm1 - i ; 00202 int ix = n3f * i ; 00203 int ixsym = n3f * isym ; 00204 // ... F+(theta) 00205 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ; 00206 // ... F_(theta) sin(theta) 00207 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ; 00208 g.set(i) = fp + fms ; 00209 g.set(isym) = fp - fms ; 00210 } 00211 //... cas particuliers: 00212 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00213 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00214 00215 // Developpement de G en series de Fourier par une FFT 00216 //---------------------------------------------------- 00217 00218 fftw_execute(p) ; 00219 00220 // Coefficients pairs du developmt. sin(l theta) de f 00221 //---------------------------------------------------- 00222 // Ces coefficients sont egaux aux coefficients en sinus du developpement 00223 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation 00224 // de fftw) : 00225 00226 double fac = -2. / double(nm1) ; 00227 cf0[0] = 0. ; 00228 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ; 00229 cf0[n3c*nm1] = 0. ; 00230 00231 // Coefficients impairs du developmt. en sin(l theta) de f 00232 //--------------------------------------------------------- 00233 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00234 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00235 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00236 // remplacer par un -2.) 00237 00238 cf0[n3c] = -fac * g(0); 00239 fac *= -2. ; 00240 for ( i = 3; i < nt; i += 2 ) { 00241 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ; 00242 } 00243 00244 } // fin de la boucle sur r 00245 00246 //------------------------------------------------------------------------ 00247 // partie sin(m phi) avec m pair : transformation en sin(l theta) 00248 //------------------------------------------------------------------------ 00249 00250 j++ ; 00251 00252 if ( j != n1f-1 ) { 00253 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont 00254 // pas nuls 00255 00256 for (k=0; k<n3f; k++) { 00257 00258 int i0 = n2n3f * j + k ; // indice de depart 00259 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00260 00261 i0 = n2n3c * j + k ; // indice de depart 00262 double* cf0 = cf + i0 ; // tableau resultat 00263 00264 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta) 00265 //--------------------------------------------- 00266 for ( i = 1; i < nm1s2 ; i++ ) { 00267 int isym = nm1 - i ; 00268 int ix = n3f * i ; 00269 int ixsym = n3f * isym ; 00270 // ... F+(theta) 00271 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ; 00272 // ... F_(theta) sin(theta) 00273 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ; 00274 g.set(i) = fp + fms ; 00275 g.set(isym) = fp - fms ; 00276 } 00277 //... cas particuliers: 00278 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00279 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00280 00281 // Developpement de G en series de Fourier par une FFT 00282 //---------------------------------------------------- 00283 00284 fftw_execute(p) ; 00285 00286 // Coefficients pairs du developmt. sin(l theta) de f 00287 //---------------------------------------------------- 00288 // Ces coefficients sont egaux aux coefficients en sinus du developpement 00289 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation 00290 // de fftw) : 00291 00292 double fac = -2. / double(nm1) ; 00293 cf0[0] = 0. ; 00294 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ; 00295 cf0[n3c*nm1] = 0. ; 00296 00297 // Coefficients impairs du developmt. en sin(l theta) de f 00298 //--------------------------------------------------------- 00299 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00300 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00301 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00302 // remplacer par un -2.) 00303 00304 cf0[n3c] = -fac * g(0); 00305 fac *= -2. ; 00306 for ( i = 3; i < nt; i += 2 ) { 00307 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ; 00308 } 00309 00310 } // fin de la boucle sur r 00311 00312 } // fin du cas ou le calcul etait necessaire (i.e. ou les 00313 // coef en phi n'etaient pas nuls) 00314 00315 00316 // On passe au cas m impair suivant: 00317 j+=3 ; 00318 00319 } // fin de la boucle sur les cas m pair 00320 00321 if (n1f<=3) // cas m=0 seulement (symetrie axiale) 00322 return ; 00323 00324 //======================================================================= 00325 // Cas m impair 00326 //======================================================================= 00327 00328 j = 2 ; 00329 00330 while (j<n1f-1) { //le dernier coef en phi n'est pas considere 00331 // (car nul) 00332 00333 //------------------------------------------------------------------------ 00334 // partie cos(m phi) avec m impair : transformation en cos(l theta) 00335 //------------------------------------------------------------------------ 00336 00337 00338 for (k=0; k<n3f; k++) { 00339 00340 int i0 = n2n3f * j + k ; // indice de depart 00341 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00342 00343 i0 = n2n3c * j + k ; // indice de depart 00344 double* cf0 = cf + i0 ; // tableau resultat 00345 00346 00347 // Valeur en theta=0 de la partie antisymetrique de F, F_ : 00348 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ); 00349 00350 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 00351 //--------------------------------------------- 00352 for ( i = 1; i < nm1s2 ; i++ ) { 00353 int isym = nm1 - i ; 00354 int ix = n3f * i ; 00355 int ixsym = n3f * isym ; 00356 // ... F+(theta) 00357 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ; 00358 // ... F_(theta) sin(psi) 00359 double fms = 0.5 * ( ff0[ix] - ff0[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 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00365 g.set(nm1s2) = ff0[ n3f*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. cos(l theta) de f 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) : 00377 00378 double fac = 2./double(nm1) ; 00379 cf0[0] = g(0) / double(nm1) ; 00380 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ; 00381 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ; 00382 00383 // Coefficients impairs du developmt. en cos(l theta) de f 00384 //--------------------------------------------------------- 00385 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00386 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00387 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00388 // remplacer par un -2.) 00389 fac *= 2. ; 00390 cf0[n3c] = 0 ; 00391 double som = 0; 00392 for ( i = 3; i < nt; i += 2 ) { 00393 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ; 00394 som += cf0[n3c*i] ; 00395 } 00396 00397 // 2. Calcul de c_1 : 00398 double c1 = ( fmoins0 - som ) / nm1s2 ; 00399 00400 // 3. Coef. c_k avec k impair: 00401 cf0[n3c] = c1 ; 00402 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ; 00403 00404 00405 } // fin de la boucle sur r 00406 00407 //-------------------------------------------------------------------- 00408 // partie sin(m phi) avec m impair : transformation en cos(l theta) 00409 //-------------------------------------------------------------------- 00410 00411 j++ ; 00412 00413 if ( (j != 1) && (j != n1f-1 ) ) { 00414 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont 00415 // pas nuls 00416 00417 for (k=0; k<n3f; k++) { 00418 00419 int i0 = n2n3f * j + k ; // indice de depart 00420 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00421 00422 i0 = n2n3c * j + k ; // indice de depart 00423 double* cf0 = cf + i0 ; // tableau resultat 00424 00425 00426 // Valeur en theta=0 de la partie antisymetrique de F, F_ : 00427 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ); 00428 00429 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 00430 //--------------------------------------------- 00431 for ( i = 1; i < nm1s2 ; i++ ) { 00432 int isym = nm1 - i ; 00433 int ix = n3f * i ; 00434 int ixsym = n3f * isym ; 00435 // ... F+(theta) 00436 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ; 00437 // ... F_(theta) sin(psi) 00438 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 00439 g.set(i) = fp + fms ; 00440 g.set(isym) = fp - fms ; 00441 } 00442 //... cas particuliers: 00443 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00444 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00445 00446 // Developpement de G en series de Fourier par une FFT 00447 //---------------------------------------------------- 00448 00449 fftw_execute(p) ; 00450 00451 // Coefficients pairs du developmt. cos(l theta) de f 00452 //---------------------------------------------------- 00453 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00454 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00455 // de fftw) : 00456 00457 double fac = 2./double(nm1) ; 00458 cf0[0] = g(0) / double(nm1) ; 00459 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ; 00460 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ; 00461 00462 // Coefficients impairs du developmt. en cos(l theta) de f 00463 //--------------------------------------------------------- 00464 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00465 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00466 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00467 // remplacer par un -2.) 00468 fac *= 2. ; 00469 cf0[n3c] = 0 ; 00470 double som = 0; 00471 for ( i = 3; i < nt; i += 2 ) { 00472 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ; 00473 som += cf0[n3c*i] ; 00474 } 00475 00476 // 2. Calcul de c_1 : 00477 double c1 = ( fmoins0 - som ) / nm1s2 ; 00478 00479 // 3. Coef. c_k avec k impair: 00480 cf0[n3c] = c1 ; 00481 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ; 00482 00483 00484 } // fin de la boucle sur r 00485 } // fin du cas ou le calcul etait necessaire (i.e. ou les 00486 // coef en phi n'etaient pas nuls) 00487 00488 // On passe au cas m impair suivant: 00489 j+=3 ; 00490 00491 } // fin de la boucle sur les cas m pair 00492 00493 }
1.4.6