cfrchebp.C

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

Generated on Tue Feb 7 01:35:15 2012 for LORENE by  doxygen 1.4.6