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