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 cftcosi_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcosi.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $" ; 00024 00025 /* 00026 * Transformation en cos((2*l+1)*theta) sur le deuxieme indice (theta) 00027 * d'un tableau 3-D representant une fonction antisymetrique par rapport 00028 * au plan z=0. 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 theta est nt = deg[1] et doit etre de la forme 00036 * nt = 2*p + 1 00037 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois 00038 * dimensions. 00039 * On doit avoir dimf[1] >= deg[1] = nt. 00040 * NB: pour dimf[0] = 1 (un seul point en phi), la transformation 00041 * est bien effectuee. 00042 * pour dimf[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 != dimf[0]-1 (cf. commentaires sur borne_phi). 00045 * 00046 * double* ff : tableau des valeurs de la fonction aux nt points de 00047 * de collocation 00048 * 00049 * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1 00050 * 00051 * L'espace memoire correspondant a ce 00052 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit 00053 * etre alloue avant l'appel a la routine. 00054 * Les valeurs de la fonction doivent etre stokees 00055 * dans le tableau ff comme suit 00056 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ] 00057 * ou j et k sont les indices correspondant a 00058 * phi et r respectivement. 00059 * 00060 * int* dimc : tableau du nombre d'elements de cc dans chacune des trois 00061 * dimensions. 00062 * On doit avoir dimc[1] >= deg[1] = nt. 00063 * Sortie: 00064 * ------- 00065 * double* cf : tableau des coefficients c_l de la fonction definis 00066 * comme suit (a r et phi fixes) 00067 * 00068 * f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) . 00069 * 00070 * L'espace memoire correspondant a ce 00071 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit 00072 * etre alloue avant l'appel a la routine. 00073 * Le coefficient c_l (0 <= l <= nt-1) est stoke dans 00074 * le tableau cf comme suit 00075 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ] 00076 * ou j et k sont les indices correspondant a 00077 * phi et r respectivement. On a c_{nt-1}=0. 00078 * 00079 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un 00080 * seul tableau, qui constitue une entree/sortie. 00081 * 00082 */ 00083 00084 /* 00085 * $Id: cftcosi.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00086 * $Log: cftcosi.C,v $ 00087 * Revision 1.1 2004/12/21 17:06:02 j_novak 00088 * Added all files for using fftw3. 00089 * 00090 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon 00091 * Suppressed the directive #include <malloc.h> for malloc is defined 00092 * in <stdlib.h> 00093 * 00094 * Revision 1.3 2002/10/16 14:36:44 j_novak 00095 * Reorganization of #include instructions of standard C++, in order to 00096 * use experimental version 3 of gcc. 00097 * 00098 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon 00099 * Modification of declaration of Fortran 77 prototypes for 00100 * a better portability (in particular on IBM AIX systems): 00101 * All Fortran subroutine names are now written F77_* and are 00102 * defined in the new file C++/Include/proto_f77.h. 00103 * 00104 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon 00105 * LORENE 00106 * 00107 * Revision 2.0 1999/02/22 15:47:57 hyc 00108 * *** empty log message *** 00109 * 00110 * 00111 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcosi.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00112 * 00113 */ 00114 00115 00116 // headers du C 00117 #include <stdlib.h> 00118 #include <fftw3.h> 00119 00120 //Lorene prototypes 00121 #include "tbl.h" 00122 00123 // Prototypage des sous-routines utilisees: 00124 fftw_plan prepare_fft(int, Tbl*&) ; 00125 double* cheb_ini(const int) ; 00126 double* chebimp_ini(const int ) ; 00127 //***************************************************************************** 00128 00129 void cftcosi(const int* deg, const int* dimf, double* ff, const int* dimc, 00130 double* cf) 00131 { 00132 00133 int i, j, k ; 00134 00135 // Dimensions des tableaux ff et cf : 00136 int n1f = dimf[0] ; 00137 int n2f = dimf[1] ; 00138 int n3f = dimf[2] ; 00139 int n1c = dimc[0] ; 00140 int n2c = dimc[1] ; 00141 int n3c = dimc[2] ; 00142 00143 // Nombre de degres de liberte en theta : 00144 int nt = deg[1] ; 00145 00146 // Tests de dimension: 00147 if (nt > n2f) { 00148 cout << "cftcosi: nt > n2f : nt = " << nt << " , n2f = " 00149 << n2f << endl ; 00150 abort () ; 00151 exit(-1) ; 00152 } 00153 if (nt > n2c) { 00154 cout << "cftcosi: nt > n2c : nt = " << nt << " , n2c = " 00155 << n2c << endl ; 00156 abort () ; 00157 exit(-1) ; 00158 } 00159 if (n1f > n1c) { 00160 cout << "cftcosi: n1f > n1c : n1f = " << n1f << " , n1c = " 00161 << n1c << endl ; 00162 abort () ; 00163 exit(-1) ; 00164 } 00165 if (n3f > n3c) { 00166 cout << "cftcosi: n3f > n3c : n3f = " << n3f << " , n3c = " 00167 << n3c << endl ; 00168 abort () ; 00169 exit(-1) ; 00170 } 00171 00172 // Nombre de points pour la FFT: 00173 int nm1 = nt - 1; 00174 int nm1s2 = nm1 / 2; 00175 00176 // Recherche des tables pour la FFT: 00177 Tbl* pg = 0x0 ; 00178 fftw_plan p = prepare_fft(nm1, pg) ; 00179 Tbl& g = *pg ; 00180 00181 // Recherche de la table des sin(psi) : 00182 double* sinp = cheb_ini(nt); 00183 00184 // Recherche de la table des points de collocations x_k : 00185 double* x = chebimp_ini(nt); 00186 00187 // boucle sur phi et r (les boucles vont resp. de 0 a max(dimf[0]-2,0) et 00188 // 0 a dimf[2]-1 ) 00189 00190 int n2n3f = n2f * n3f ; 00191 int n2n3c = n2c * n3c ; 00192 00193 /* 00194 * Borne de la boucle sur phi: 00195 * si n1f = 1, on effectue la boucle une fois seulement. 00196 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients 00197 * j=n1f-1 et j=0 ne sont pas consideres car nuls). 00198 */ 00199 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ; 00200 00201 for (j=0; j< borne_phi; j++) { 00202 00203 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi) 00204 00205 for (k=0; k<n3f; k++) { 00206 00207 int i0 = n2n3f * j + k ; // indice de depart 00208 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00209 00210 i0 = n2n3c * j + k ; // indice de depart 00211 double* cf0 = cf + i0 ; // tableau resultat 00212 00213 // Multiplication de la fonction par x=cos(theta) (pour la rendre paire) 00214 // NB: dans les commentaires qui suivent, on note h(x) = x f(x). 00215 // (Les valeurs de h dans l'espace des configurations sont stokees dans le 00216 // tableau cf0). 00217 for (i=0; i<nt-1; i++) cf0[n3c*i] = x[nm1-i] * ff0[n3f*i] ; 00218 cf0[n3c*nm1] = 0 ; 00219 00220 /* 00221 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi] 00222 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)). 00223 */ 00224 00225 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00226 double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] ); 00227 00228 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 00229 //--------------------------------------------- 00230 for ( i = 1; i < nm1s2 ; i++ ) { 00231 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2: 00232 int isym = nm1 - i ; 00233 // ... indice (dans le tableau cf0) du point theta correspondant a psi 00234 int ix = n3c * i ; 00235 // ... indice (dans le tableau cf0) du point theta correspondant a sym(psi) 00236 int ixsym = n3c * isym ; 00237 // ... F+(psi) 00238 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ; 00239 // ... F_(psi) sin(psi) 00240 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ; 00241 g.set(i) = fp + fms ; 00242 g.set(isym) = fp - fms ; 00243 } 00244 //... cas particuliers: 00245 g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] ); 00246 g.set(nm1s2) = cf0[ n3c*nm1s2 ]; 00247 00248 // Developpement de G en series de Fourier par une FFT 00249 //---------------------------------------------------- 00250 00251 fftw_execute(p) ; 00252 00253 // Coefficients pairs du developmt. de Tchebyshev de h 00254 //---------------------------------------------------- 00255 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00256 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00257 // de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le 00258 // remplacer par un +1.) : 00259 00260 double fac = 2./double(nm1) ; 00261 cf0[0] = g(0) / double(nm1) ; 00262 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ; 00263 cf0[n3c*nm1] = g(nm1s2) ; 00264 00265 // Coefficients impairs du developmt. de Tchebyshev de h 00266 //------------------------------------------------------ 00267 // 1. Coef. c'_k (recurrence amorcee a partir de zero) 00268 // Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw 00269 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00270 // remplacer par un -2.) 00271 fac *= 2. ; 00272 cf0[n3c] = 0 ; 00273 double som = 0; 00274 for ( i = 3; i < nt; i += 2 ) { 00275 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ; 00276 som += cf0[n3c*i] ; 00277 } 00278 00279 // 2. Calcul de c_1 : 00280 double c1 = ( fmoins0 - som ) / nm1s2 ; 00281 00282 // 3. Coef. c_k avec k impair: 00283 cf0[n3c] = c1 ; 00284 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ; 00285 00286 00287 // Coefficients de f en fonction de ceux de h 00288 //------------------------------------------- 00289 00290 cf0[0] = 2* cf0[0] ; 00291 for (i=1; i<nm1; i++) { 00292 cf0[n3c*i] = 2 * cf0[n3c*i] - cf0[n3c*(i-1)] ; 00293 } 00294 cf0[n3c*nm1] = 0 ; 00295 00296 } // fin de la boucle sur r 00297 } // fin de la boucle sur phi 00298 00299 }
1.4.6