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

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