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