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 cftcossinc_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinc.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 cos(l theta ) . 00069 * pour m impair: 00070 * f(theta) = som_{l=0}^{nt-1} c_l sin( 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: cftcossinc.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00090 * $Log: cftcossinc.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/cftcossinc.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 00107 //Lorene prototypes 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 cftcossinc(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 << "cftcossinc: nt > n2f : nt = " << nt << " , n2f = " 00135 << n2f << endl ; 00136 abort () ; 00137 exit(-1) ; 00138 } 00139 if (nt > n2c) { 00140 cout << "cftcossinc: nt > n2c : nt = " << nt << " , n2c = " 00141 << n2c << endl ; 00142 abort () ; 00143 exit(-1) ; 00144 } 00145 if (n1f > n1c) { 00146 cout << "cftcossinc: n1f > n1c : n1f = " << n1f << " , n1c = " 00147 << n1c << endl ; 00148 abort () ; 00149 exit(-1) ; 00150 } 00151 if (n3f > n3c) { 00152 cout << "cftcossinc: 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 cos(l theta) 00187 //------------------------------------------------------------------------ 00188 00189 00190 for (k=0; k<n3f; k++) { 00191 00192 int i0 = n2n3f * j + k ; // indice de depart 00193 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00194 00195 i0 = n2n3c * j + k ; // indice de depart 00196 double* cf0 = cf + i0 ; // tableau resultat 00197 00198 00199 // Valeur en theta=0 de la partie antisymetrique de F, F_ : 00200 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ); 00201 00202 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 00203 //--------------------------------------------- 00204 for ( i = 1; i < nm1s2 ; i++ ) { 00205 int isym = nm1 - i ; 00206 int ix = n3f * i ; 00207 int ixsym = n3f * isym ; 00208 // ... F+(theta) 00209 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ; 00210 // ... F_(theta) sin(psi) 00211 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 00212 g.set(i) = fp + fms ; 00213 g.set(isym) = fp - fms ; 00214 } 00215 //... cas particuliers: 00216 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00217 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00218 00219 // Developpement de G en series de Fourier par une FFT 00220 //---------------------------------------------------- 00221 00222 fftw_execute(p) ; 00223 00224 // Coefficients pairs du developmt. cos(l theta) de f 00225 //---------------------------------------------------- 00226 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00227 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00228 // de fftw) : 00229 00230 double fac = 2./double(nm1) ; 00231 cf0[0] = g(0) / double(nm1) ; 00232 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ; 00233 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ; 00234 00235 // Coefficients impairs du developmt. en cos(l theta) de f 00236 //--------------------------------------------------------- 00237 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00238 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00239 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00240 // remplacer par un -2.) 00241 fac *= 2. ; 00242 cf0[n3c] = 0 ; 00243 double som = 0; 00244 for ( i = 3; i < nt; i += 2 ) { 00245 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ; 00246 som += cf0[n3c*i] ; 00247 } 00248 00249 // 2. Calcul de c_1 : 00250 double c1 = ( fmoins0 - som ) / nm1s2 ; 00251 00252 // 3. Coef. c_k avec k impair: 00253 cf0[n3c] = c1 ; 00254 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ; 00255 00256 00257 } // fin de la boucle sur r 00258 00259 //-------------------------------------------------------------------- 00260 // partie sin(m phi) avec m pair : transformation en cos(l theta) 00261 //-------------------------------------------------------------------- 00262 00263 j++ ; 00264 00265 if ( (j != 1) && (j != n1f-1 ) ) { 00266 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont 00267 // pas nuls 00268 00269 for (k=0; k<n3f; k++) { 00270 00271 int i0 = n2n3f * j + k ; // indice de depart 00272 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00273 00274 i0 = n2n3c * j + k ; // indice de depart 00275 double* cf0 = cf + i0 ; // tableau resultat 00276 00277 00278 // Valeur en theta=0 de la partie antisymetrique de F, F_ : 00279 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ); 00280 00281 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 00282 //--------------------------------------------- 00283 for ( i = 1; i < nm1s2 ; i++ ) { 00284 int isym = nm1 - i ; 00285 int ix = n3f * i ; 00286 int ixsym = n3f * isym ; 00287 // ... F+(theta) 00288 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ; 00289 // ... F_(theta) sin(psi) 00290 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 00291 g.set(i) = fp + fms ; 00292 g.set(isym) = fp - fms ; 00293 } 00294 //... cas particuliers: 00295 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00296 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00297 00298 // Developpement de G en series de Fourier par une FFT 00299 //---------------------------------------------------- 00300 00301 fftw_execute(p) ; 00302 00303 // Coefficients pairs du developmt. cos(l theta) de f 00304 //---------------------------------------------------- 00305 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00306 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00307 // de fftw) : 00308 00309 double fac = 2./double(nm1) ; 00310 cf0[0] = g(0) / double(nm1) ; 00311 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ; 00312 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ; 00313 00314 // Coefficients impairs du developmt. en cos(l theta) de f 00315 //--------------------------------------------------------- 00316 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00317 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00318 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00319 // remplacer par un -2.) 00320 fac *= 2. ; 00321 cf0[n3c] = 0 ; 00322 double som = 0; 00323 for ( i = 3; i < nt; i += 2 ) { 00324 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ; 00325 som += cf0[n3c*i] ; 00326 } 00327 00328 // 2. Calcul de c_1 : 00329 double c1 = ( fmoins0 - som ) / nm1s2 ; 00330 00331 // 3. Coef. c_k avec k impair: 00332 cf0[n3c] = c1 ; 00333 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ; 00334 00335 00336 } // fin de la boucle sur r 00337 } // fin du cas ou le calcul etait necessaire (i.e. ou les 00338 // coef en phi n'etaient pas nuls) 00339 00340 // On passe au cas m pair suivant: 00341 j+=3 ; 00342 00343 } // fin de la boucle sur les cas m pair 00344 00345 if (n1f<=3) // cas m=0 seulement (symetrie axiale) 00346 return ; 00347 00348 //======================================================================= 00349 // Cas m impair 00350 //======================================================================= 00351 00352 j = 2 ; 00353 00354 while (j<n1f-1) { //le dernier coef en phi n'est pas considere 00355 // (car nul) 00356 00357 //-------------------------------------------------------------------- 00358 // partie cos(m phi) avec m impair : transformation en sin(l) theta) 00359 //-------------------------------------------------------------------- 00360 00361 for (k=0; k<n3f; k++) { 00362 00363 int i0 = n2n3f * j + k ; // indice de depart 00364 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00365 00366 i0 = n2n3c * j + k ; // indice de depart 00367 double* cf0 = cf + i0 ; // tableau resultat 00368 00369 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta) 00370 //--------------------------------------------- 00371 00372 for ( i = 1; i < nm1s2 ; i++ ) { 00373 int isym = nm1 - i ; 00374 int ix = n3f * i ; 00375 int ixsym = n3f * isym ; 00376 // ... F+(theta) 00377 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ; 00378 // ... F_(theta) sin(theta) 00379 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ; 00380 g.set(i) = fp + fms ; 00381 g.set(isym) = fp - fms ; 00382 } 00383 //... cas particuliers: 00384 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00385 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00386 00387 // Developpement de G en series de Fourier par une FFT 00388 //---------------------------------------------------- 00389 00390 fftw_execute(p) ; 00391 00392 // Coefficients pairs du developmt. sin(l theta) de f 00393 //---------------------------------------------------- 00394 // Ces coefficients sont egaux aux coefficients en sinus du developpement 00395 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation 00396 // de fftw) : 00397 00398 double fac = -2. / double(nm1) ; 00399 cf0[0] = 0. ; 00400 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ; 00401 cf0[n3c*nm1] = 0. ; 00402 00403 // Coefficients impairs du developmt. en sin(l theta) de f 00404 //--------------------------------------------------------- 00405 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00406 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00407 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00408 // remplacer par un -2.) 00409 00410 cf0[n3c] = -fac * g(0); 00411 fac *= -2. ; 00412 for ( i = 3; i < nt; i += 2 ) { 00413 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ; 00414 } 00415 00416 } // fin de la boucle sur r 00417 00418 //------------------------------------------------------------------------ 00419 // partie sin(m phi) avec m impair : transformation en sin(l theta) 00420 //------------------------------------------------------------------------ 00421 00422 j++ ; 00423 00424 if ( j != n1f-1 ) { 00425 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont 00426 // pas nuls 00427 00428 for (k=0; k<n3f; k++) { 00429 00430 int i0 = n2n3f * j + k ; // indice de depart 00431 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00432 00433 i0 = n2n3c * j + k ; // indice de depart 00434 double* cf0 = cf + i0 ; // tableau resultat 00435 00436 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta) 00437 //--------------------------------------------- 00438 for ( i = 1; i < nm1s2 ; i++ ) { 00439 int isym = nm1 - i ; 00440 int ix = n3f * i ; 00441 int ixsym = n3f * isym ; 00442 // ... F+(theta) 00443 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ; 00444 // ... F_(theta) sin(theta) 00445 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ; 00446 g.set(i) = fp + fms ; 00447 g.set(isym) = fp - fms ; 00448 } 00449 //... cas particuliers: 00450 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] ); 00451 g.set(nm1s2) = ff0[ n3f*nm1s2 ]; 00452 00453 // Developpement de G en series de Fourier par une FFT 00454 //---------------------------------------------------- 00455 00456 fftw_execute(p) ; 00457 00458 // Coefficients pairs du developmt. sin(l theta) de f 00459 //---------------------------------------------------- 00460 // Ces coefficients sont egaux aux coefficients en sinus du developpement 00461 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation 00462 // de fftw) : 00463 00464 double fac = -2. / double(nm1) ; 00465 cf0[0] = 0. ; 00466 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ; 00467 cf0[n3c*nm1] = 0. ; 00468 00469 // Coefficients impairs du developmt. en sin(l theta) de f 00470 //--------------------------------------------------------- 00471 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00472 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw 00473 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00474 // remplacer par un -2.) 00475 00476 cf0[n3c] = -fac * g(0); 00477 fac *= -2. ; 00478 for ( i = 3; i < nt; i += 2 ) { 00479 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ; 00480 } 00481 00482 } // fin de la boucle sur r 00483 00484 } // fin du cas ou le calcul etait necessaire (i.e. ou les 00485 // coef en phi n'etaient pas nuls) 00486 00487 00488 // On passe au cas m impair suivant: 00489 j+=3 ; 00490 00491 } // fin de la boucle sur les cas m impair 00492 00493 }
1.4.6