cftcossincp.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 cftcossincp_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossincp.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $" ;
00024 
00025 
00026 /*
00027  * Transformation en cos(2*l*theta) ou sin((2*l+1)*theta) (suivant la parite
00028  *  de l'indice m en phi) sur le deuxieme indice (theta)
00029  *  d'un tableau 3-D representant une fonction symetrique par rapport
00030  *  au plan z=0.
00031  *  Utilise la bibliotheque fftw.
00032  *
00033  * Entree:
00034  * -------
00035  *   int* deg   : tableau du nombre effectif de degres de liberte dans chacune
00036  *        des 3 dimensions: le nombre de points de collocation
00037  *        en theta est  nt = deg[1] et doit etre de la forme
00038  *          nt = 2*p + 1 
00039  *   int* dimf  : tableau du nombre d'elements de ff dans chacune des trois 
00040  *            dimensions.
00041  *        On doit avoir  dimf[1] >= deg[1] = nt. 
00042  *
00043  *   double* ff : tableau des valeurs de la fonction aux nt points de
00044  *                        de collocation
00045  *
00046  *            theta_l =  pi/2 l/(nt-1)       0 <= l <= nt-1 
00047  *
00048  *            L'espace memoire correspondant a ce
00049  *                        pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit 
00050  *            etre alloue avant l'appel a la routine.    
00051  *            Les valeurs de la fonction doivent etre stokees
00052  *            dans le tableau ff comme suit
00053  *          f( theta_l ) = ff[ dimf[1]*dimf[2] * m + k + dimf[2] * l ]
00054  *           ou m et k sont les indices correspondant a
00055  *           phi et r respectivement.
00056  *  NB: cette routine suppose que la transformation en phi a deja ete
00057  *      effectuee: ainsi m est un indice de Fourier, non un indice de
00058  *      point de collocation en phi.
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  *            pour m pair:
00069  *          f(theta) = som_{l=0}^{nt-1} c_l cos( 2 l theta ) . 
00070  *            pour m impair:
00071  *          f(theta) = som_{l=0}^{nt-2} c_l sin( (2 l+1) theta ) . 
00072  *
00073  *            L'espace memoire correspondant a ce
00074  *                        pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit 
00075  *            etre alloue avant l'appel a la routine.    
00076  *           Le coefficient c_l (0 <= l <= nt-1) est stoke dans
00077  *           le tableau cf comme suit
00078  *                c_l = cf[ dimc[1]*dimc[2] * m + k + dimc[2] * l ]
00079  *           ou m et k sont les indices correspondant a
00080  *           phi et r respectivement.
00081  *           Pour m impair, c_{nt-1} = 0.
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 /*
00089  * $Id: cftcossincp.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $
00090  * $Log: cftcossincp.C,v $
00091  * Revision 1.1  2004/12/21 17:06:02  j_novak
00092  * Added all files for using fftw3.
00093  *
00094  * Revision 1.4  2003/01/31 10:31:23  e_gourgoulhon
00095  * Suppressed the directive #include <malloc.h> for malloc is defined
00096  * in <stdlib.h>
00097  *
00098  * Revision 1.3  2002/10/16 14:36:51  j_novak
00099  * Reorganization of #include instructions of standard C++, in order to
00100  * use experimental version 3 of gcc.
00101  *
00102  * Revision 1.2  2002/09/09 13:00:39  e_gourgoulhon
00103  * Modification of declaration of Fortran 77 prototypes for
00104  * a better portability (in particular on IBM AIX systems):
00105  * All Fortran subroutine names are now written F77_* and are
00106  * defined in the new file C++/Include/proto_f77.h.
00107  *
00108  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00109  * LORENE
00110  *
00111  * Revision 2.1  2000/01/27  12:16:02  eric
00112  * Modif commentaires.
00113  *
00114  * Revision 2.0  1999/02/22  15:47:32  hyc
00115  * *** empty log message ***
00116  *
00117  *
00118  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossincp.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $
00119  *
00120  */
00121 
00122 
00123 // headers du C
00124 #include <stdlib.h>
00125 #include <fftw3.h>
00126 
00127 //Lorene prototypes
00128 #include "tbl.h"
00129 
00130 // Prototypage des sous-routines utilisees:
00131 fftw_plan prepare_fft(int, Tbl*&) ;
00132 double* cheb_ini(const int) ;
00133 double* chebimp_ini(const int ) ;
00134 //*****************************************************************************
00135 
00136 void cftcossincp(const int* deg, const int* dimf, double* ff, const int* dimc,
00137         double* cf)
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 // Nombre de degres de liberte en theta :    
00151     int nt = deg[1] ;
00152     
00153 // Tests de dimension:
00154     if (nt > n2f) {
00155     cout << "cftcossincp: nt > n2f : nt = " << nt << " ,  n2f = " 
00156     << n2f << endl ;
00157     abort () ;
00158     exit(-1) ;
00159     }
00160     if (nt > n2c) {
00161     cout << "cftcossincp: nt > n2c : nt = " << nt << " ,  n2c = " 
00162     << n2c << endl ;
00163     abort () ;
00164     exit(-1) ;
00165     }
00166     if (n1f > n1c) {
00167     cout << "cftcossincp: n1f > n1c : n1f = " << n1f << " ,  n1c = " 
00168     << n1c << endl ;
00169     abort () ;
00170     exit(-1) ;
00171     }
00172     if (n3f > n3c) {
00173     cout << "cftcossincp: n3f > n3c : n3f = " << n3f << " ,  n3c = " 
00174     << n3c << endl ;
00175     abort () ;
00176     exit(-1) ;
00177     }
00178 
00179 // Nombre de points pour la FFT:
00180     int nm1 = nt - 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(nt);    
00190     
00191 // Recherche de la table des sin( theta_l ) :
00192     double* sinth = chebimp_ini(nt);    
00193 
00194 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
00195 //           et 0 a dimf[2])
00196 
00197     int n2n3f = n2f * n3f ;
00198     int n2n3c = n2c * n3c ;
00199 
00200 //=======================================================================
00201 //              Cas m pair
00202 //=======================================================================
00203 
00204     j = 0 ;
00205     
00206     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00207             // (car nul)
00208 
00209 //--------------------------------------------------------------------
00210 //  partie cos(m phi) avec m pair : transformation en cos(2 l theta)
00211 //--------------------------------------------------------------------
00212 
00213     for (k=0; k<n3f; k++) {
00214 
00215         int i0 = n2n3f * j + k ; // indice de depart 
00216         double* ff0 = ff + i0 ;    // tableau des donnees a transformer
00217 
00218         i0 = n2n3c * j + k ; // indice de depart 
00219         double* cf0 = cf + i0 ;    // tableau resultat
00220 
00221 /*
00222  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00223  *     reliee a theta par  psi = 2 theta   et F(psi) = f(theta(psi)).  
00224  */
00225  
00226 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00227             double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
00228 
00229 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 
00230 //---------------------------------------------
00231             for ( i = 1; i < nm1s2 ; i++ ) {
00232 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00233         int isym = nm1 - i ; 
00234 // ... indice (dans le tableau ff0) du point theta correspondant a psi
00235         int ix = n3f * i ;
00236 // ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
00237         int ixsym = n3f * isym ;
00238 // ... F+(psi)
00239         double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;    
00240 // ... F_(psi) sin(psi)
00241         double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 
00242         g.set(i) = fp + fms ;
00243         g.set(isym) = fp - fms ;
00244             }
00245 //... cas particuliers:
00246             g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
00247             g.set(nm1s2) = ff0[ n3f*nm1s2 ];
00248 
00249 // Developpement de G en series de Fourier par une FFT
00250 //----------------------------------------------------
00251 
00252         fftw_execute(p) ;
00253 
00254 // Coefficients pairs du developmt. cos(2l theta) de f
00255 //----------------------------------------------------
00256 //  Ces coefficients sont egaux aux coefficients en cosinus du developpement
00257 //  de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
00258 //  de fftw) :
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) / double(nm1) ;    
00264 
00265 // Coefficients impairs du developmt. en cos(2l theta) de f
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     }   // fin de la boucle sur r 
00288 
00289 //--------------------------------------------------------------------
00290 //  partie sin(m phi) avec m pair : transformation en cos(2 l theta)
00291 //--------------------------------------------------------------------
00292 
00293     j++ ;
00294 
00295     if ( (j != 1) && (j != n1f-1 ) ) {  
00296 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00297 //  pas nuls 
00298     for (k=0; k<n3f; k++) {
00299 
00300         int i0 = n2n3f * j + k ; // indice de depart 
00301         double* ff0 = ff + i0 ;    // tableau des donnees a transformer
00302 
00303         i0 = n2n3c * j + k ; // indice de depart 
00304         double* cf0 = cf + i0 ;    // tableau resultat
00305 
00306 /*
00307  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00308  *     reliee a theta par  psi = 2 theta   et F(psi) = f(theta(psi)).  
00309  */
00310  
00311 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00312             double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
00313 
00314 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 
00315 //---------------------------------------------
00316             for ( i = 1; i < nm1s2 ; i++ ) {
00317 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00318         int isym = nm1 - i ; 
00319 // ... indice (dans le tableau ff0) du point theta correspondant a psi
00320         int ix = n3f * i ;
00321 // ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
00322         int ixsym = n3f * isym ;
00323 // ... F+(psi)
00324         double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;    
00325 // ... F_(psi) sin(psi)
00326         double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ; 
00327         g.set(i) = fp + fms ;
00328         g.set(isym) = fp - fms ;
00329             }
00330 //... cas particuliers:
00331             g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
00332             g.set(nm1s2) = ff0[ n3f*nm1s2 ];
00333 
00334 // Developpement de G en series de Fourier par une FFT
00335 //----------------------------------------------------
00336 
00337         fftw_execute(p) ;
00338 
00339 // Coefficients pairs du developmt. cos(2l theta) de f
00340 //----------------------------------------------------
00341 //  Ces coefficients sont egaux aux coefficients en cosinus du developpement
00342 //  de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
00343 //  de fftw) :
00344 
00345         double fac = 2./double(nm1) ;
00346         cf0[0] = g(0) / double(nm1) ;
00347             for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ; 
00348         cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;    
00349 
00350 // Coefficients impairs du developmt. en cos(2l theta) de f
00351 //---------------------------------------------------------
00352 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
00353 //  Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
00354 //  (si fftw donnait reellement les coef. en sinus, il faudrait le
00355 //   remplacer par un -2.) 
00356         fac *= 2. ;
00357             cf0[n3c] = 0 ;
00358             double som = 0;
00359             for ( i = 3; i < nt; i += 2 ) {
00360         cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
00361             som += cf0[n3c*i] ;
00362             }
00363 
00364 // 2. Calcul de c_1 :
00365             double c1 = ( fmoins0 - som ) / nm1s2 ;
00366 
00367 // 3. Coef. c_k avec k impair:  
00368             cf0[n3c] = c1 ;
00369             for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
00370 
00371 
00372     }   // fin de la boucle sur r 
00373     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00374          // coef en phi n'etaient pas nuls)
00375 
00376 // On passe au cas m pair suivant:
00377     j+=3 ;
00378 
00379    }    // fin de la boucle sur les cas m pair
00380 
00381     if (n1f<=3)         // cas m=0 seulement (symetrie axiale)
00382     return ;
00383     
00384 //=======================================================================
00385 //              Cas m impair
00386 //=======================================================================
00387 
00388     j = 2 ;
00389     
00390     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00391             // (car nul)
00392 
00393 //------------------------------------------------------------------------
00394 //  partie cos(m phi) avec m impair : transformation en sin((2 l+1) theta)
00395 //------------------------------------------------------------------------
00396 
00397     for (k=0; k<n3f; k++) {
00398 
00399         int i0 = n2n3f * j + k ; // indice de depart 
00400         double* ff0 = ff + i0 ;    // tableau des donnees a transformer
00401 
00402         i0 = n2n3c * j + k ; // indice de depart 
00403         double* cf0 = cf + i0 ;    // tableau resultat
00404 
00405 // Multiplication de la fonction par sin(theta) (pour la rendre developpable
00406 //  en cos(2l theta) ) 
00407 // NB: dans les commentaires qui suivent, on note 
00408 //     h(theta) = f(theta) sin(theta).
00409 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
00410 //  tableau cf0).
00411         cf0[0] = 0 ;
00412         for (i=1; i<nt; i++) cf0[n3c*i] = sinth[i] * ff0[n3f*i] ;
00413 
00414 /*
00415  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00416  *     reliee a theta par  psi = 2 theta  et F(psi) = h(theta(psi)).  
00417  */
00418  
00419 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00420             double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
00421 
00422 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 
00423 //---------------------------------------------
00424             for ( i = 1; i < nm1s2 ; i++ ) {
00425 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00426         int isym = nm1 - i ; 
00427 // ... indice (dans le tableau cf0) du point theta correspondant a psi
00428         int ix = n3c * i ;
00429 // ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
00430         int ixsym = n3c * isym ;
00431 // ... F+(psi)
00432         double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;    
00433 // ... F_(psi) sin(psi)
00434         double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ; 
00435         g.set(i) = fp + fms ;
00436         g.set(isym) = fp - fms ;
00437             }
00438 //... cas particuliers:
00439             g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
00440             g.set(nm1s2) = cf0[ n3c*nm1s2 ];
00441 
00442 // Developpement de G en series de Fourier par une FFT
00443 //----------------------------------------------------
00444 
00445         fftw_execute(p) ;
00446 
00447 // Coefficients pairs du developmt. cos(2l theta) de h
00448 //----------------------------------------------------
00449 //  Ces coefficients sont egaux aux coefficients en cosinus du developpement
00450 //  de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
00451 //  de fftw) :
00452 
00453         double fac = 2./double(nm1) ;
00454         cf0[0] = g(0)/double(nm1) ;
00455             for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;   
00456         cf0[n3c*nm1] = g(nm1s2)/double(nm1) ;    
00457 
00458 // Coefficients impairs du developmt. en cos(2l theta) de h
00459 //---------------------------------------------------------
00460 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
00461 //  Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
00462 //  (si fftw donnait reellement les coef. en sinus, il faudrait le
00463 //   remplacer par un -2.) 
00464         fac *= 2. ;
00465             cf0[n3c] = 0 ;
00466             double som = 0;
00467             for ( i = 3; i < nt; i += 2 ) {
00468         cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
00469             som += cf0[n3c*i] ;
00470             }
00471 
00472 // 2. Calcul de c_1 :
00473             double c1 = ( fmoins0 - som ) / nm1s2 ;
00474 
00475 // 3. Coef. c_k avec k impair:  
00476             cf0[n3c] = c1 ;
00477             for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
00478 
00479 // Coefficients de f en fonction de ceux de h
00480 //-------------------------------------------
00481 
00482             cf0[0] = 2* cf0[0] ;
00483             for (i=1; i<nm1; i++) {
00484         cf0[n3c*i] = 2 * cf0[n3c*i] + cf0[n3c*(i-1)] ;
00485             }    
00486             cf0[n3c*nm1] = 0 ;
00487 
00488     }   // fin de la boucle sur r 
00489 
00490 //------------------------------------------------------------------------
00491 //  partie sin(m phi) avec m impair : transformation en sin((2 l+1) theta)
00492 //------------------------------------------------------------------------
00493 
00494     j++ ;
00495 
00496     if ( j != n1f-1  ) {  
00497 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00498 //  pas nuls 
00499 
00500     for (k=0; k<n3f; k++) {
00501 
00502         int i0 = n2n3f * j + k ; // indice de depart 
00503         double* ff0 = ff + i0 ;    // tableau des donnees a transformer
00504 
00505         i0 = n2n3c * j + k ; // indice de depart 
00506         double* cf0 = cf + i0 ;    // tableau resultat
00507 
00508 // Multiplication de la fonction par sin(theta) (pour la rendre developpable
00509 //  en cos(2l theta) ) 
00510 // NB: dans les commentaires qui suivent, on note 
00511 //     h(theta) = f(theta) sin(theta).
00512 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
00513 //  tableau cf0).
00514         cf0[0] = 0 ;
00515         for (i=1; i<nt; i++) cf0[n3c*i] = sinth[i] * ff0[n3f*i] ;
00516 
00517 /*
00518  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00519  *     reliee a theta par  psi = 2 theta  et F(psi) = h(theta(psi)).  
00520  */
00521  
00522 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00523             double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
00524 
00525 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 
00526 //---------------------------------------------
00527             for ( i = 1; i < nm1s2 ; i++ ) {
00528 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00529         int isym = nm1 - i ; 
00530 // ... indice (dans le tableau cf0) du point theta correspondant a psi
00531         int ix = n3c * i ;
00532 // ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
00533         int ixsym = n3c * isym ;
00534 // ... F+(psi)
00535         double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;    
00536 // ... F_(psi) sin(psi)
00537         double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ; 
00538         g.set(i) = fp + fms ;
00539         g.set(isym) = fp - fms ;
00540             }
00541 //... cas particuliers:
00542             g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
00543             g.set(nm1s2) = cf0[ n3c*nm1s2 ];
00544 
00545 // Developpement de G en series de Fourier par une FFT
00546 //----------------------------------------------------
00547 
00548         fftw_execute(p) ;
00549 
00550 // Coefficients pairs du developmt. cos(2l theta) de h
00551 //----------------------------------------------------
00552 //  Ces coefficients sont egaux aux coefficients en cosinus du developpement
00553 //  de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
00554 //  de fftw) :
00555 
00556         double fac = 2./double(nm1) ;
00557         cf0[0] = g(0)/double(nm1) ;
00558             for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;   
00559         cf0[n3c*nm1] = g(nm1s2)/double(nm1) ;    
00560 
00561 // Coefficients impairs du developmt. en cos(2l theta) de h
00562 //---------------------------------------------------------
00563 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
00564 //  Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
00565 //  (si fftw donnait reellement les coef. en sinus, il faudrait le
00566 //   remplacer par un -2.) 
00567         fac *= 2. ;
00568             cf0[n3c] = 0 ;
00569             double som = 0;
00570             for ( i = 3; i < nt; i += 2 ) {
00571         cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
00572             som += cf0[n3c*i] ;
00573             }
00574 
00575 // 2. Calcul de c_1 :
00576             double c1 = ( fmoins0 - som ) / nm1s2 ;
00577 
00578 // 3. Coef. c_k avec k impair:  
00579             cf0[n3c] = c1 ;
00580             for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
00581 
00582 // Coefficients de f en fonction de ceux de h
00583 //-------------------------------------------
00584 
00585             cf0[0] = 2* cf0[0] ;
00586             for (i=1; i<nm1; i++) {
00587         cf0[n3c*i] = 2 * cf0[n3c*i] + cf0[n3c*(i-1)] ;
00588             }    
00589             cf0[n3c*nm1] = 0 ;
00590 
00591     }   // fin de la boucle sur r 
00592 
00593     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00594          // coef en phi n'etaient pas nuls)
00595 
00596 
00597 // On passe au cas m impair suivant:
00598     j+=3 ;
00599 
00600    }    // fin de la boucle sur les cas m impair
00601 
00602 }

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