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 cfrcheb_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfrcheb.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $" ; 00024 00025 /* 00026 * Transformation de Tchebyshev (cas fin) sur le troisieme indice (indice 00027 * correspondant a r) d'un tableau 3-D 00028 * par le biais de la bibliotheque fftw. 00029 * 00030 * Entree: 00031 * ------- 00032 * int* deg : tableau du nombre effectif de degres de liberte dans chacune 00033 * des 3 dimensions: le nombre de points de collocation 00034 * en r est nr = deg[2] et doit etre de la forme 00035 * nr = 2*p + 1 00036 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois 00037 * dimensions. 00038 * On doit avoir dimf[2] >= deg[2] = nr. 00039 * NB: pour dimf[0] = 1 (un seul point en phi), la transformation 00040 * est bien effectuee. 00041 * pour dimf[0] > 1 (plus d'un point en phi), la 00042 * transformation n'est effectuee que pour les indices (en phi) 00043 * j != 1 et j != dimf[0]-1 (cf. commentaires sur borne_phi). 00044 * 00045 * double* ff : tableau des valeurs de la fonction aux nr points de 00046 * de collocation 00047 * 00048 * x_i = - cos( pi i/(nr-1) ) 0 <= i <= nr-1 00049 * 00050 * Les valeurs de la fonction doivent etre stokees dans le 00051 * tableau ff comme suit 00052 * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ] 00053 * ou j et k sont les indices correspondant a phi et theta 00054 * respectivement. 00055 * L'espace memoire correspondant a ce pointeur doit etre 00056 * dimf[0]*dimf[1]*dimf[2] et doit avoir ete alloue avant 00057 * l'appel a la routine. 00058 * 00059 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois 00060 * dimensions. 00061 * On doit avoir dimc[2] >= deg[2] = nr. 00062 * 00063 * Sortie: 00064 * ------- 00065 * double* cf : tableau des coefficients c_i de la fonction definis 00066 * comme suit (a theta et phi fixes) 00067 * 00068 * f(x) = som_{i=0}^{nr-1} c_i T_i(x) , 00069 * 00070 * ou T_i(x) designe le polynome de Tchebyshev de degre i. 00071 * Les coefficients c_i (0 <= i <= nr-1) sont stokes dans 00072 * le tableau cf comme suit 00073 * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ] 00074 * ou j et k sont les indices correspondant a phi et theta 00075 * respectivement. 00076 * L'espace memoire correspondant au pointeur cf doit etre 00077 * dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant 00078 * l'appel a la routine. 00079 * 00080 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un 00081 * seul tableau, qui constitue une entree/sortie. 00082 */ 00083 00084 /* 00085 * $Id: cfrcheb.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00086 * $Log: cfrcheb.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:43 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:48:47 hyc 00108 * *** empty log message *** 00109 * 00110 * 00111 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfrcheb.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 <assert.h> 00119 #include <fftw3.h> 00120 00121 //Lorene prototypes 00122 #include "tbl.h" 00123 00124 // Prototypage des sous-routines utilisees: 00125 fftw_plan prepare_fft(int, Tbl*&) ; 00126 double* cheb_ini(const int) ; 00127 //***************************************************************************** 00128 00129 void cfrcheb(const int* deg, const int* dimf, double* ff, const int* dimc, 00130 double* cf) 00131 00132 { 00133 00134 int i, j, k ; 00135 00136 // Dimensions des tableaux ff et cf : 00137 int n1f = dimf[0] ; 00138 int n2f = dimf[1] ; 00139 int n3f = dimf[2] ; 00140 int n1c = dimc[0] ; 00141 int n2c = dimc[1] ; 00142 int n3c = dimc[2] ; 00143 00144 // Nombres de degres de liberte en r : 00145 int nr = deg[2] ; 00146 00147 // Tests de dimension: 00148 if (nr > n3f) { 00149 cout << "cfrcheb: nr > n3f : nr = " << nr << " , n3f = " 00150 << n3f << endl ; 00151 abort () ; 00152 exit(-1) ; 00153 } 00154 if (nr > n3c) { 00155 cout << "cfrcheb: nr > n3c : nr = " << nr << " , n3c = " 00156 << n3c << endl ; 00157 abort () ; 00158 exit(-1) ; 00159 } 00160 if (n1f > n1c) { 00161 cout << "cfrcheb: n1f > n1c : n1f = " << n1f << " , n1c = " 00162 << n1c << endl ; 00163 abort () ; 00164 exit(-1) ; 00165 } 00166 if (n2f > n2c) { 00167 cout << "cfrcheb: n2f > n2c : n2f = " << n2f << " , n2c = " 00168 << n2c << endl ; 00169 abort () ; 00170 exit(-1) ; 00171 } 00172 00173 00174 // Nombre de points pour la FFT: 00175 int nm1 = nr - 1; 00176 int nm1s2 = nm1 / 2; 00177 00178 // Recherche des tables pour la FFT: 00179 Tbl* pg = 0x0 ; 00180 fftw_plan p = prepare_fft(nm1, pg) ; 00181 Tbl& g = *pg ; 00182 00183 // Recherche de la table des sin(psi) : 00184 double* sinp = cheb_ini(nr); 00185 00186 // boucle sur phi et theta 00187 00188 int n2n3f = n2f * n3f ; 00189 int n2n3c = n2c * n3c ; 00190 00191 /* 00192 * Borne de la boucle sur phi: 00193 * si n1f = 1, on effectue la boucle une fois seulement. 00194 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients 00195 * j=n1f-1 et j=0 ne sont pas consideres car nuls). 00196 */ 00197 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ; 00198 00199 for (j=0; j< borne_phi; j++) { 00200 00201 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi) 00202 00203 for (k=0; k<n2f; k++) { 00204 00205 int i0 = n2n3f * j + n3f * k ; // indice de depart 00206 double* ff0 = ff + i0 ; // tableau des donnees a transformer 00207 00208 i0 = n2n3c * j + n3c * k ; // indice de depart 00209 double* cf0 = cf + i0 ; // tableau resultat 00210 00211 /* 00212 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi] 00213 * reliee a x par x = - cos(psi) et F(psi) = f(x(psi)). 00214 */ 00215 00216 // Valeur en psi=0 de la partie antisymetrique de F, F_ : 00217 double fmoins0 = 0.5 * ( ff0[0] - ff0[nm1] ); 00218 00219 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 00220 //--------------------------------------------- 00221 for ( i = 1; i < nm1s2 ; i++ ) { 00222 // ... indice du pt symetrique de psi par rapport a pi/2: 00223 int isym = nm1 - i ; 00224 // ... F+(psi) 00225 double fp = 0.5 * ( ff0[i] + ff0[isym] ) ; 00226 // ... F_(psi) sin(psi) 00227 double fms = 0.5 * ( ff0[i] - ff0[isym] ) * sinp[i] ; 00228 g.set(i) = fp + fms ; 00229 g.set(isym) = fp - fms ; 00230 } 00231 //... cas particuliers: 00232 g.set(0) = 0.5 * ( ff0[0] + ff0[nm1] ); 00233 g.set(nm1s2) = ff0[nm1s2]; 00234 00235 // Developpement de G en series de Fourier par une FFT 00236 //---------------------------------------------------- 00237 00238 fftw_execute(p) ; 00239 00240 // Coefficients pairs du developmt. de Tchebyshev de f 00241 //---------------------------------------------------- 00242 // Ces coefficients sont egaux aux coefficients en cosinus du developpement 00243 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation 00244 // de fftw) : 00245 00246 double fac = 2./double(nm1) ; 00247 cf0[0] = g(0) / double(nm1) ; 00248 for (i=2; i<nm1; i += 2) cf0[i] = fac*g(i/2) ; 00249 cf0[nm1] = g(nm1s2) / double(nm1) ; 00250 00251 // Coefficients impairs du developmt. de Tchebyshev de f 00252 //------------------------------------------------------ 00253 // 1. Coef. c'_k (recurrence amorcee a partir de zero): 00254 // NB: Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw 00255 // (si fftw donnait reellement les coef. en sinus, il faudrait le 00256 // remplacer par un +2.) 00257 fac *= -2. ; 00258 cf0[1] = 0 ; 00259 double som = 0; 00260 for ( i = 3; i < nr; i += 2 ) { 00261 cf0[i] = cf0[i-2] + fac * g(nm1-i/2) ; 00262 som += cf0[i] ; 00263 } 00264 00265 // 2. Calcul de c_1 : 00266 double c1 = - ( fmoins0 + som ) / nm1s2 ; 00267 00268 // 3. Coef. c_k avec k impair: 00269 cf0[1] = c1 ; 00270 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ; 00271 00272 00273 } // fin de la boucle sur theta 00274 } // fin de la boucle sur phi 00275 00276 }
1.4.6