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 cfrchebpip_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfrchebpip.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $" ; 00024 00025 00026 /* 00027 * Transformation de Tchebyshev (cas rare) sur le troisieme indice (indice 00028 * correspondant a r) d'un tableau 3-D decrivant une fonction quelconque. 00029 * Utilise la bibliotheque fftw. 00030 * 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 r est nr = deg[2] et doit etre de la forme 00037 * nr = 2*p + 1 00038 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois 00039 * dimensions. 00040 * On doit avoir dimf[2] >= deg[2] = nr. 00041 * NB: pour dimf[0] = 1 (un seul point en phi), la transformation 00042 * est bien effectuee. 00043 * pour dimf[0] > 1 (plus d'un point en phi), la 00044 * transformation n'est effectuee que pour les indices (en phi) 00045 * j != 1 et j != dimf[0]-1 (cf. commentaires sur borne_phi). 00046 * 00047 * double* ff : tableau des valeurs de la fonction aux nr points de 00048 * de collocation 00049 * 00050 * x_i = sin( pi/2 i/(nr-1) ) 0 <= i <= nr-1 00051 * 00052 * Les valeurs de la fonction doivent etre stokees dans le 00053 * tableau ff comme suit 00054 * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ] 00055 * ou j et k sont les indices correspondant a phi et theta 00056 * respectivement. 00057 * L'espace memoire correspondant a ce pointeur doit etre 00058 * dimf[0]*dimf[1]*dimf[2] et doit etre alloue avant l'appel a 00059 * la routine. 00060 * 00061 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois 00062 * dimensions. 00063 * On doit avoir dimc[2] >= deg[2] = nr. 00064 * 00065 * Sortie: 00066 * ------- 00067 * double* cf : tableau des nr coefficients c_i de la fonction definis 00068 * comme suit (a theta et phi fixes) 00069 * 00070 * Si l pair f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) , 00071 * Si l impair f(x) = som_{i=0}^{nr-1} c_i T_{2i+1}(x) 00072 * 00073 * ou T_{i}(x) designe le polynome de Tchebyshev de degre i. 00074 * Les coefficients c_i (0 <= i <= nr-1) sont stokes dans 00075 * le tableau cf comme suit 00076 * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ] 00077 * ou j et k sont les indices correspondant a phi et theta 00078 * respectivement. 00079 * L'espace memoire correspondant a ce pointeur doit etre 00080 * dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant 00081 * l'appel a la routine. 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 * $Id: cfrchebpip.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00089 * $Log: cfrchebpip.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.1 2004/11/23 15:13:50 m_forot 00094 * Added the bases for the cases without any equatorial symmetry 00095 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I). 00096 * 00097 * 00098 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfrchebpip.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00099 * 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 double* chebimp_ini(const int ) ; 00114 00115 //***************************************************************************** 00116 00117 void cfrchebpip(const int* deg, const int* dimf, double* ff, const int* dimc, 00118 double* cf) 00119 00120 { 00121 00122 int i, j, k ; 00123 00124 // Dimensions des tableaux ff et cf : 00125 int n1f = dimf[0] ; 00126 int n2f = dimf[1] ; 00127 int n3f = dimf[2] ; 00128 int n1c = dimc[0] ; 00129 int n2c = dimc[1] ; 00130 int n3c = dimc[2] ; 00131 00132 // Nombres de degres de liberte en r : 00133 int nr = deg[2] ; 00134 00135 // Tests de dimension: 00136 if (nr > n3f) { 00137 cout << "cfrchebpip: nr > n3f : nr = " << nr << " , n3f = " 00138 << n3f << endl ; 00139 abort () ; 00140 exit(-1) ; 00141 } 00142 if (nr > n3c) { 00143 cout << "cfrchebpip: nr > n3c : nr = " << nr << " , n3c = " 00144 << n3c << endl ; 00145 abort () ; 00146 exit(-1) ; 00147 } 00148 if (n1f > n1c) { 00149 cout << "cfrchebpip: n1f > n1c : n1f = " << n1f << " , n1c = " 00150 << n1c << endl ; 00151 abort () ; 00152 exit(-1) ; 00153 } 00154 if (n2f > n2c) { 00155 cout << "cfrchebpip: n2f > n2c : n2f = " << n2f << " , n2c = " 00156 << n2c << endl ; 00157 abort () ; 00158 exit(-1) ; 00159 } 00160 00161 // Nombre de points pour la FFT: 00162 int nm1 = nr - 1; 00163 int nm1s2 = nm1 / 2; 00164 00165 // Recherche des tables pour la FFT: 00166 Tbl* pg = 0x0 ; 00167 fftw_plan p = prepare_fft(nm1, pg) ; 00168 Tbl& g = *pg ; 00169 00170 // Recherche de la table des sin(psi) : 00171 double* sinp = cheb_ini(nr); 00172 00173 // Recherche de la table des points de collocations x_k : 00174 double* x = chebimp_ini(nr); 00175 00176 // boucle sur phi et theta 00177 00178 int n2n3f = n2f * n3f ; 00179 int n2n3c = n2c * n3c ; 00180 00181 /* 00182 * Borne de la boucle sur phi: 00183 * si n1f = 1, on effectue la boucle une fois seulement. 00184 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients 00185 * j=n1f-1 et j=0 ne sont pas consideres car nuls). 00186 */ 00187 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ; 00188 00189 for (j=0; j< borne_phi; j++) { 00190 00191 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi) 00192 00193 /************** Cas l pair ***************/ 00194 00195 for (k=0; k<n2f; k+=2) { 00196 00197 int i0 = n2n3f * j + n3f * k ; // indice de depart 00198 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00199 00200 i0 = n2n3c * j + n3c * k ; // indice de depart 00201 double* cf0 = cf + i0 ; // tableau resultat 00202 00203 00204 // Valeur en theta=0 de la partie antisymetrique de F, F_ : 00205 double fmoins0 = 0.5 * ( ff0[nm1] - ff0[0] ); 00206 00207 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 00208 //--------------------------------------------- 00209 for ( i = 1; i < nm1s2 ; i++ ) { 00210 // ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2: 00211 int isym = nm1 - i ; 00212 // ... indice (dans le tableau ff0) du point x correspondant a theta 00213 int ix = nm1 - i ; 00214 // ... indice (dans le tableau ff0) du point x correspondant a sym(theta) 00215 int ixsym = nm1 - isym ; 00216 00217 // ... F+(psi) 00218 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ; 00219 // ... F_(psi) sin(psi) 00220 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 00221 g.set(i) = fp + fms ; 00222 g.set(isym) = fp - fms ; 00223 } 00224 //... cas particuliers: 00225 g.set(0) = 0.5 * ( ff0[nm1] + ff0[0] ); 00226 g.set(nm1s2) = ff0[nm1s2]; 00227 00228 // Developpement de G en series de Fourier par une FFT 00229 //---------------------------------------------------- 00230 00231 fftw_execute(p) ; 00232 00233 // Coefficients pairs du developmt. de Tchebyshev de f 00234 //---------------------------------------------------- 00235 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00236 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00237 // de fftw) : 00238 00239 double fac = 2./double(nm1) ; 00240 cf0[0] = g(0) / double(nm1) ; 00241 for (i=2; i<nm1; i += 2) cf0[i] = fac*g(i/2) ; 00242 cf0[nm1] = g(nm1s2) / double(nm1) ; 00243 00244 // Coefficients impairs du developmt. de Tchebyshev de f 00245 //------------------------------------------------------ 00246 // 1. Coef. c'_k (recurrence amorcee a partir de zero) 00247 // Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw 00248 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00249 // remplacer par un -2.) 00250 fac *= 2. ; 00251 cf0[1] = 0. ; 00252 double som = 0.; 00253 for ( i = 3; i < nr; i += 2 ) { 00254 cf0[i] = cf0[i-2] + fac * g(nm1-i/2) ; 00255 som += cf0[i] ; 00256 } 00257 00258 // 2. Calcul de c_1 : 00259 double c1 = ( fmoins0 - som ) / nm1s2 ; 00260 00261 // 3. Coef. c_k avec k impair: 00262 cf0[1] = c1 ; 00263 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ; 00264 00265 00266 } // fin de la boucle sur theta 00267 00268 /************** Cas l impair ***************/ 00269 00270 for (k=1; k<n2f; k+=2) { 00271 int i0 = n2n3f * j + n3f * k ; // indice de depart 00272 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00273 00274 i0 = n2n3c * j + n3c * k ; // indice de depart 00275 double* cf0 = cf + i0 ; // tableau resultat 00276 00277 // Multiplication de la fonction par x (pour la rendre paire) 00278 // NB: dans les commentaires qui suivent, on note h(x) = x f(x). 00279 // (Les valeurs de h dans l'espace des configurations sont stokees dans le 00280 // tableau cf0). 00281 cf0[0] = 0 ; 00282 for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ; 00283 00284 00285 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00286 double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] ); 00287 00288 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 00289 //--------------------------------------------- 00290 for ( i = 1; i < nm1s2 ; i++ ) { 00291 // ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2: 00292 int isym = nm1 - i ; 00293 // ... indice (dans le tableau cf0) du point x correspondant a theta 00294 int ix = nm1 - i ; 00295 // ... indice (dans le tableau cf0) du point x correspondant a sym(theta) 00296 int ixsym = nm1 - isym ; 00297 00298 // ... F+(psi) 00299 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ; 00300 // ... F_(psi) sin(psi) 00301 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ; 00302 g.set(i) = fp + fms ; 00303 g.set(isym) = fp - fms ; 00304 } 00305 //... cas particuliers: 00306 g.set(0) = 0.5 * ( cf0[nm1] + cf0[0] ); 00307 g.set(nm1s2) = cf0[nm1s2]; 00308 00309 // Developpement de G en series de Fourier par une FFT 00310 //---------------------------------------------------- 00311 00312 fftw_execute(p) ; 00313 00314 // Coefficients pairs du developmt. de Tchebyshev de h 00315 //---------------------------------------------------- 00316 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00317 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00318 // de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le 00319 // remplacer par un +1.) : 00320 00321 double fac = 2./double(nm1) ; 00322 cf0[0] = g(0) / double(nm1) ; 00323 for (i=2; i<nm1; i += 2 ) cf0[i] = fac * g(i/2) ; 00324 cf0[nm1] = g(nm1s2) / double(nm1) ; 00325 00326 // Coefficients impairs du developmt. de Tchebyshev de h 00327 //------------------------------------------------------ 00328 // 1. Coef. c'_k (recurrence amorcee a partir de zero) 00329 // Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw 00330 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00331 // remplacer par un -2.) 00332 fac *= 2 ; 00333 cf0[1] = 0 ; 00334 double som = 0; 00335 for ( i = 3; i < nr; i += 2 ) { 00336 cf0[i] = cf0[i-2] + fac * g(nm1 - i/2) ; 00337 som += cf0[i] ; 00338 } 00339 00340 // 2. Calcul de c_1 : 00341 double c1 = ( fmoins0 - som ) / nm1s2 ; 00342 00343 // 3. Coef. c_k avec k impair: 00344 cf0[1] = c1 ; 00345 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ; 00346 00347 // Coefficients de f en fonction de ceux de h 00348 //------------------------------------------- 00349 00350 cf0[0] = 2* cf0[0] ; 00351 for (i=1; i<nm1; i++) { 00352 cf0[i] = 2 * cf0[i] - cf0[i-1] ; 00353 } 00354 cf0[nm1] = 0 ; 00355 00356 00357 } // fin de la boucle sur theta 00358 } // fin de la boucle sur phi 00359 00360 }
1.4.6