circhebpimi.C

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

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