cfrchebpii.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 cfrchebpii_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cfrchebpii.C,v 1.1 2004/12/21 17:06:01 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 quelconque. 
00029  *  Utilise la routine FFT Fortran FFT991
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 3^q 5^r + 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  *           Si l impair  f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) ,
00071  *                   Si l pair    f(x) = som_{i=0}^{nr-1} c_i T_{2i+1}(x) 
00072  *
00073  *          ou T_{i}(x) designe le polynome de Tchebyshev de degre i.    
00074  *          Les coefficients c_i (0 <= i <= nr-1) sont stokes dans
00075  *          le tableau cf comme suit
00076  *             c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ]
00077  *          ou j et k sont les indices correspondant a phi et theta 
00078  *          respectivement.
00079  *          L'espace memoire correspondant a ce pointeur doit etre 
00080  *          dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant 
00081  *          l'appel a la routine.    
00082  *
00083  * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un 
00084  *     seul tableau, qui constitue une entree/sortie.
00085  */
00086  
00087 /*
00088  * $Id: cfrchebpii.C,v 1.1 2004/12/21 17:06:01 j_novak Exp $
00089  * $Log: cfrchebpii.C,v $
00090  * Revision 1.1  2004/12/21 17:06:01  j_novak
00091  * Added all files for using fftw3.
00092  *
00093  * Revision 1.1  2004/11/23 15:13:50  m_forot
00094  * Added the bases for the cases without any equatorial symmetry
00095  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
00096  *
00097  *
00098  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cfrchebpii.C,v 1.1 2004/12/21 17:06:01 j_novak Exp $
00099  *
00100  */
00101 
00102 
00103 // headers du C
00104 #include <assert.h>
00105 #include <stdlib.h>
00106 
00107 // Prototypes of F77 subroutines
00108 #include "headcpp.h"
00109 #include "proto_f77.h"
00110 
00111 // Prototypage des sous-routines utilisees:
00112 int*    facto_ini(int ) ;
00113 double* trigo_ini(int ) ;
00114 double* cheb_ini(const int) ;
00115 double* chebimp_ini(const int ) ;
00116 
00117 //*****************************************************************************
00118 
00119 void cfrchebpii(const int* deg, const int* dimf, double* ff, const int* dimc,
00120         double* cf)
00121 
00122 {
00123 
00124 int i, j, k ;
00125 
00126 // Dimensions des tableaux ff et cf  :
00127     int n1f = dimf[0] ;
00128     int n2f = dimf[1] ;
00129     int n3f = dimf[2] ;
00130     int n1c = dimc[0] ;
00131     int n2c = dimc[1] ;
00132     int n3c = dimc[2] ;
00133 
00134 // Nombres de degres de liberte en r :    
00135     int nr = deg[2] ;
00136     
00137 // Tests de dimension:
00138     if (nr > n3f) {
00139     cout << "cfrchebpii: nr > n3f : nr = " << nr << " ,  n3f = " 
00140     << n3f << endl ;
00141     abort () ;
00142     exit(-1) ;
00143     }
00144     if (nr > n3c) {
00145     cout << "cfrchebpii: nr > n3c : nr = " << nr << " ,  n3c = " 
00146     << n3c << endl ;
00147     abort () ;
00148     exit(-1) ;
00149     }
00150     if (n1f > n1c) {
00151     cout << "cfrchebpii: n1f > n1c : n1f = " << n1f << " ,  n1c = " 
00152     << n1c << endl ;
00153     abort () ;
00154     exit(-1) ;
00155     }
00156     if (n2f > n2c) {
00157     cout << "cfrchebpii: n2f > n2c : n2f = " << n2f << " ,  n2c = " 
00158     << n2c << endl ;
00159     abort () ;
00160     exit(-1) ;
00161     }
00162 
00163 // Nombre de points pour la FFT:
00164     int nm1 = nr - 1;
00165     int nm1s2 = nm1 / 2;
00166 
00167 // Recherche des tables pour la FFT:
00168     int* facto = facto_ini(nm1) ;
00169     double* trigo = trigo_ini(nm1) ;
00170 
00171 // Recherche de la table des sin(psi) :
00172     double* sinp = cheb_ini(nr);    
00173 
00174 // Recherche de la table des points de collocations x_k :
00175     double* x = chebimp_ini(nr);
00176     
00177     // tableau de travail G et t1
00178     // (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
00179     double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) );  
00180     double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
00181 
00182 // Parametres pour la routine FFT991
00183     int jump = 1 ;
00184     int inc = 1 ;
00185     int lot = 1 ;
00186     int isign = - 1 ;
00187 
00188 // boucle sur phi et theta
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     /************** Cas l impair ***************/
00206 
00207     for (k=1; k<n2f; k+=2) {
00208 
00209         int i0 = n2n3f * j + n3f * k ; // indice de depart 
00210         double* ff0 = ff + i0 ;    // tableau des donnees a transformer
00211 
00212         i0 = n2n3c * j + n3c * k ; // indice de depart 
00213         double* cf0 = cf + i0 ;    // tableau resultat
00214 
00215  
00216 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
00217             double fmoins0 = 0.5 * ( ff0[nm1] - ff0[0] );
00218 
00219 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 
00220 //---------------------------------------------
00221             for ( i = 1; i < nm1s2 ; i++ ) {
00222 // ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
00223         int isym = nm1 - i ; 
00224 // ... indice (dans le tableau ff0) du point x correspondant a theta
00225         int ix = nm1 - i ;
00226 // ... indice (dans le tableau ff0) du point x correspondant a sym(theta)
00227         int ixsym = nm1 -  isym ;
00228     
00229 // ... F+(psi)
00230         double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;    
00231 // ... F_(psi) sin(psi)
00232         double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 
00233         g[i] = fp + fms ;
00234         g[isym] = fp - fms ;
00235             }
00236 //... cas particuliers:
00237             g[0] = 0.5 * ( ff0[nm1] + ff0[0] );
00238             g[nm1s2] = ff0[nm1s2];
00239 
00240 // Developpement de G en series de Fourier par une FFT
00241 //----------------------------------------------------
00242 
00243             F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
00244 
00245 // Coefficients pairs du developmt. de Tchebyshev de f
00246 //----------------------------------------------------
00247 //  Ces coefficients sont egaux aux coefficients en cosinus du developpement
00248 //  de G en series de Fourier (le facteur 2 vient de la normalisation
00249 //  de fft991) :
00250 
00251         cf0[0] = g[0] ;
00252             for (i=2; i<nm1; i += 2 ) cf0[i] = 2. * g[i] ;  
00253             cf0[nm1] = g[nm1] ;
00254 
00255 // Coefficients impairs du developmt. de Tchebyshev de f
00256 //------------------------------------------------------
00257 // 1. Coef. c'_k (recurrence amorcee a partir de zero)
00258 //  Le +4. en facteur de g[i] est du a la normalisation de fft991
00259 //  (si fft991 donnait reellement les coef. en sinus, il faudrait le
00260 //   remplacer par un -2.) 
00261             cf0[1] = 0 ;
00262             double som = 0;
00263             for ( i = 3; i < nr; i += 2 ) {
00264         cf0[i] = cf0[i-2] + 4. * g[i] ;
00265             som += cf0[i] ;
00266             }
00267 
00268 // 2. Calcul de c_1 :
00269             double c1 = ( fmoins0 - som ) / nm1s2 ;
00270 
00271 // 3. Coef. c_k avec k impair:  
00272             cf0[1] = c1 ;
00273             for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
00274 
00275 
00276     }   // fin de la boucle sur theta
00277 
00278     /************** Cas l pair ***************/
00279 
00280     for (k=0; k<n2f; k+=2) {
00281         int i0 = n2n3f * j + n3f * k ; // indice de depart 
00282         double* ff0 = ff + i0 ;    // tableau des donnees a transformer
00283 
00284         i0 = n2n3c * j + n3c * k ; // indice de depart 
00285         double* cf0 = cf + i0 ;    // tableau resultat
00286 
00287 // Multiplication de la fonction par x (pour la rendre paire) 
00288 // NB: dans les commentaires qui suivent, on note h(x) = x f(x).
00289 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
00290 //  tableau cf0).
00291         cf0[0] = 0 ;
00292         for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
00293 
00294 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00295             double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
00296 
00297 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta) 
00298 //---------------------------------------------
00299             for ( i = 1; i < nm1s2 ; i++ ) {
00300 // ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
00301         int isym = nm1 - i ; 
00302 // ... indice (dans le tableau cf0) du point x correspondant a theta
00303         int ix = nm1 - i ;
00304 // ... indice (dans le tableau cf0) du point x correspondant a sym(theta)
00305         int ixsym = nm1 -  isym ;
00306     
00307 // ... F+(psi)
00308         double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;    
00309 // ... F_(psi) sin(psi)
00310         double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ; 
00311         g[i] = fp + fms ;
00312         g[isym] = fp - fms ;
00313             }
00314 //... cas particuliers:
00315             g[0] = 0.5 * ( cf0[nm1] + cf0[0] );
00316             g[nm1s2] = cf0[nm1s2];
00317 
00318 // Developpement de G en series de Fourier par une FFT
00319 //----------------------------------------------------
00320 
00321             F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
00322 
00323 // Coefficients pairs du developmt. de Tchebyshev de h
00324 //----------------------------------------------------
00325 //  Ces coefficients sont egaux aux coefficients en cosinus du developpement
00326 //  de G en series de Fourier (le facteur 2. vient de la normalisation
00327 //  de fft991; si fft991 donnait reellement les coef. en cosinus, il faudrait le
00328 //   remplacer par un +1.)  :
00329 
00330         cf0[0] = g[0] ;
00331             for (i=2; i<nm1; i += 2 ) cf0[i] = 2. * g[i] ;  
00332             cf0[nm1] = g[nm1] ;
00333 
00334 // Coefficients impairs du developmt. de Tchebyshev de h
00335 //------------------------------------------------------
00336 // 1. Coef. c'_k (recurrence amorcee a partir de zero)
00337 //  Le +4. en facteur de g[i] est du a la normalisation de fft991
00338 //  (si fft991 donnait reellement les coef. en sinus, il faudrait le
00339 //   remplacer par un -2.) 
00340             cf0[1] = 0 ;
00341             double som = 0;
00342             for ( i = 3; i < nr; i += 2 ) {
00343         cf0[i] = cf0[i-2] + 4. * g[i] ;
00344             som += cf0[i] ;
00345             }
00346 
00347 // 2. Calcul de c_1 :
00348             double c1 = ( fmoins0 - som ) / nm1s2 ;
00349 
00350 // 3. Coef. c_k avec k impair:  
00351             cf0[1] = c1 ;
00352             for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
00353 
00354 // Coefficients de f en fonction de ceux de h
00355 //-------------------------------------------
00356 
00357     cf0[0] = 2* cf0[0] ;
00358     for (i=1; i<nm1; i++) {
00359     cf0[i] = 2 * cf0[i] - cf0[i-1] ;
00360     }    
00361     cf0[nm1] = 0 ;
00362 
00363 
00364     }   // fin de la boucle sur theta  
00365    }    // fin de la boucle sur phi
00366 
00367     // Menage
00368     free (t1) ;
00369     free (g) ;
00370 
00371 }

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