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/FFT991/circhebpimi.C,v 1.1 2004/12/21 17:06:01 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 routine FFT Fortran FFT991
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 3^q 5^r + 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.1 2004/12/21 17:06:01 j_novak Exp $
00100  * $Log: circhebpimi.C,v $
00101  * Revision 1.1  2004/12/21 17:06:01  j_novak
00102  * Added all files for using fftw3.
00103  *
00104  * Revision 1.4  2003/01/31 10:31:23  e_gourgoulhon
00105  * Suppressed the directive #include <malloc.h> for malloc is defined
00106  * in <stdlib.h>
00107  *
00108  * Revision 1.3  2002/10/16 14:36:53  j_novak
00109  * Reorganization of #include instructions of standard C++, in order to
00110  * use experimental version 3 of gcc.
00111  *
00112  * Revision 1.2  2002/09/09 13:00:40  e_gourgoulhon
00113  * Modification of declaration of Fortran 77 prototypes for
00114  * a better portability (in particular on IBM AIX systems):
00115  * All Fortran subroutine names are now written F77_* and are
00116  * defined in the new file C++/Include/proto_f77.h.
00117  *
00118  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00119  * LORENE
00120  *
00121  * Revision 2.0  1999/02/22  15:43:19  hyc
00122  * *** empty log message ***
00123  *
00124  *
00125  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/circhebpimi.C,v 1.1 2004/12/21 17:06:01 j_novak Exp $
00126  *
00127  */
00128 
00129 
00130 
00131 // headers du C
00132 #include <assert.h>
00133 #include <stdlib.h>
00134 
00135 // Prototypes of F77 subroutines
00136 #include "headcpp.h"
00137 #include "proto_f77.h"
00138 
00139 
00140 // Prototypage des sous-routines utilisees:
00141 int*    facto_ini(int ) ;
00142 double* trigo_ini(int ) ;
00143 double* cheb_ini(const int) ;
00144 double* chebimp_ini(const int ) ;
00145 //*****************************************************************************
00146 
00147 void circhebpimi(const int* deg, const int* dimc, double* cf,
00148             const int* dimf, double* ff)
00149 
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     int* facto = facto_ini(nm1) ;
00197     double* trigo = trigo_ini(nm1) ;
00198 
00199 // Recherche de la table des sin(psi) :
00200     double* sinp = cheb_ini(nr);    
00201     
00202 // Recherche de la table des points de collocations x_k :
00203     double* x = chebimp_ini(nr);    
00204 
00205     // tableau de travail t1 et g
00206     //   (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
00207     double* g =  (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;    
00208     double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
00209 
00210 // Parametres pour la routine FFT991
00211     int jump = 1 ;
00212     int inc = 1 ;
00213     int lot = 1 ;
00214     int isign = 1 ;
00215 
00216 // boucle sur phi et theta
00217 
00218     int n2n3f = n2f * n3f ;
00219     int n2n3c = n2c * n3c ;
00220 
00221 //=======================================================================
00222 //              Cas m pair
00223 //=======================================================================
00224 
00225     j = 0 ;
00226     
00227     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00228             // (car nul)
00229 
00230 
00231 //------------------------------------------------------------------------
00232 //  partie cos(m phi) avec m pair : developpement en T_{2i+1}(x)
00233 //------------------------------------------------------------------------
00234 
00235     for (k=0; k<n2c; k++) {
00236 
00237         int i0 = n2n3c * j + n3c * k ; // indice de depart 
00238         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00239 
00240         i0 = n2n3f * j + n3f * k ; // indice de depart 
00241         double* ff0 = ff + i0 ;    // tableau resultat
00242 
00243 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
00244 //  h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le
00245 //  tableau t1 :
00246     t1[0] = .5 * cf0[0] ;
00247     for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ;
00248     t1[nm1] = .5 * cf0[nr-2] ;
00249 
00250 /*
00251  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00252  *     reliee a x par  x = cos(psi/2)   et F(psi) = h(x(psi)).  
00253  */
00254 
00255 // Calcul des coefficients de Fourier de la fonction 
00256 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00257 // en fonction des coefficients de Tchebyshev de f:
00258 
00259 // Coefficients impairs de G
00260 //--------------------------
00261  
00262         double c1 = t1[1] ;
00263 
00264             double som = 0;
00265         ff0[1] = 0 ;
00266             for ( i = 3; i < nr; i += 2 ) {
00267             ff0[i] = t1[i] - c1 ;
00268         som += ff0[i] ;
00269             }   
00270 
00271 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00272         double fmoins0 = nm1s2 * c1 + som ;
00273 
00274 // Coef. impairs de G
00275 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
00276 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00277             g[1] = 0 ;
00278             for ( i = 3; i < nr; i += 2 ) {
00279         g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
00280             }
00281             g[nr] = 0 ; 
00282 
00283 
00284 // Coefficients pairs de G
00285 //------------------------
00286 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00287 //   f.
00288 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
00289 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00290 
00291         g[0] = t1[0] ;
00292             for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * t1[i] ;  
00293             g[nm1] = t1[nm1] ;
00294 
00295 // Transformation de Fourier inverse de G 
00296 //---------------------------------------
00297 
00298 // FFT inverse
00299             F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
00300 
00301 // Valeurs de f deduites de celles de G
00302 //-------------------------------------
00303 
00304             for ( i = 1; i < nm1s2 ; i++ ) {
00305 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00306         int isym = nm1 - i ; 
00307 // ... indice (dans le tableau ff0) du point x correspondant a psi
00308         int ix = nm1 - i ;
00309 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
00310         int ixsym = nm1 -  isym ;
00311 
00312         double fp = .5 * ( g[i] + g[isym] ) ;
00313         double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
00314 
00315         ff0[ix] = ( fp + fm ) / x[ix];
00316         ff0[ixsym] = ( fp - fm ) / x[ixsym] ;
00317             }
00318     
00319 //... cas particuliers:
00320         ff0[0] = 0 ;
00321         ff0[nm1] = g[0] + fmoins0 ;
00322         ff0[nm1s2] = g[nm1s2] / x[nm1s2] ;
00323 
00324     }   // fin de la boucle sur theta 
00325 
00326 //------------------------------------------------------------------------
00327 //  partie sin(m phi) avec m pair : developpement en T_{2i+1}(x)
00328 //------------------------------------------------------------------------
00329 
00330     j++ ;
00331 
00332     if ( (j != 1) && (j != n1f-1) ) {  
00333 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00334 //  pas nuls 
00335 
00336     for (k=0; k<n2c; k++) {
00337 
00338         int i0 = n2n3c * j + n3c * k ; // indice de depart 
00339         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00340 
00341         i0 = n2n3f * j + n3f * k ; // indice de depart 
00342         double* ff0 = ff + i0 ;    // tableau resultat
00343 
00344 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
00345 //  h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le
00346 //  tableau t1 :
00347     t1[0] = .5 * cf0[0] ;
00348     for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ;
00349     t1[nm1] = .5 * cf0[nr-2] ;
00350 
00351 /*
00352  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00353  *     reliee a x par  x = cos(psi/2)   et F(psi) = h(x(psi)).  
00354  */
00355 
00356 // Calcul des coefficients de Fourier de la fonction 
00357 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00358 // en fonction des coefficients de Tchebyshev de f:
00359 
00360 // Coefficients impairs de G
00361 //--------------------------
00362  
00363         double c1 = t1[1] ;
00364 
00365             double som = 0;
00366         ff0[1] = 0 ;
00367             for ( i = 3; i < nr; i += 2 ) {
00368             ff0[i] = t1[i] - c1 ;
00369         som += ff0[i] ;
00370             }   
00371 
00372 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00373         double fmoins0 = nm1s2 * c1 + som ;
00374 
00375 // Coef. impairs de G
00376 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
00377 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00378             g[1] = 0 ;
00379             for ( i = 3; i < nr; i += 2 ) {
00380         g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
00381             }
00382             g[nr] = 0 ; 
00383 
00384 
00385 // Coefficients pairs de G
00386 //------------------------
00387 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00388 //   f.
00389 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
00390 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00391 
00392         g[0] = t1[0] ;
00393             for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * t1[i] ;  
00394             g[nm1] = t1[nm1] ;
00395 
00396 // Transformation de Fourier inverse de G 
00397 //---------------------------------------
00398 
00399 // FFT inverse
00400             F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
00401 
00402 // Valeurs de f deduites de celles de G
00403 //-------------------------------------
00404 
00405             for ( i = 1; i < nm1s2 ; i++ ) {
00406 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00407         int isym = nm1 - i ; 
00408 // ... indice (dans le tableau ff0) du point x correspondant a psi
00409         int ix = nm1 - i ;
00410 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
00411         int ixsym = nm1 -  isym ;
00412 
00413         double fp = .5 * ( g[i] + g[isym] ) ;
00414         double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
00415 
00416         ff0[ix] = ( fp + fm ) / x[ix];
00417         ff0[ixsym] = ( fp - fm ) / x[ixsym] ;
00418             }
00419     
00420 //... cas particuliers:
00421         ff0[0] = 0 ;
00422         ff0[nm1] = g[0] + fmoins0 ;
00423         ff0[nm1s2] = g[nm1s2] / x[nm1s2] ;
00424 
00425     }   // fin de la boucle sur theta 
00426 
00427     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00428          // coef en phi n'etaient pas nuls)
00429 
00430 // On passe au cas m pair suivant:
00431     j+=3 ;
00432 
00433    }    // fin de la boucle sur les cas m pair
00434 
00435     if (n1f<=3) {       // cas m=0 seulement (symetrie axiale)
00436     free (g) ;
00437     free (t1) ;
00438     return ;
00439     }
00440 
00441 //=======================================================================
00442 //              Cas m impair
00443 //=======================================================================
00444 
00445     j = 2 ;
00446     
00447     while (j<n1f-1) {   
00448     
00449 //--------------------------------------------------------------------
00450 //  partie cos(m phi) avec m impair : developpement en T_{2i}(x)
00451 //--------------------------------------------------------------------
00452 
00453     for (k=0; k<n2c; k++) {
00454 
00455         int i0 = n2n3c * j + n3c * k ; // indice de depart 
00456         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00457 
00458         i0 = n2n3f * j + n3f * k ; // indice de depart 
00459         double* ff0 = ff + i0 ;    // tableau resultat
00460 
00461 /*
00462  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00463  *     reliee a x par  x = cos(psi/2)   et F(psi) = f(x(psi)).  
00464  */
00465 
00466 // Calcul des coefficients de Fourier de la fonction 
00467 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00468 // en fonction des coefficients de Tchebyshev de f:
00469 
00470 // Coefficients impairs de G
00471 //--------------------------
00472  
00473         double c1 = cf0[1] ;
00474 
00475             double som = 0;
00476         ff0[1] = 0 ;
00477             for ( i = 3; i < nr; i += 2 ) {
00478             ff0[i] = cf0[i] - c1 ;
00479         som += ff0[i] ;
00480             }   
00481 
00482 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00483         double fmoins0 = nm1s2 * c1 + som ;
00484 
00485 // Coef. impairs de G
00486 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
00487 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00488             g[1] = 0 ;
00489             for ( i = 3; i < nr; i += 2 ) {
00490         g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
00491             }
00492             g[nr] = 0 ; 
00493 
00494 
00495 // Coefficients pairs de G
00496 //------------------------
00497 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00498 //   f.
00499 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
00500 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00501 
00502         g[0] = cf0[0] ;
00503             for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * cf0[i] ; 
00504             g[nm1] = cf0[nm1] ;
00505 
00506 // Transformation de Fourier inverse de G 
00507 //---------------------------------------
00508 
00509 // FFT inverse
00510             F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
00511 
00512 // Valeurs de f deduites de celles de G
00513 //-------------------------------------
00514 
00515             for ( i = 1; i < nm1s2 ; i++ ) {
00516 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00517         int isym = nm1 - i ; 
00518 // ... indice (dans le tableau ff0) du point x correspondant a psi
00519         int ix = nm1 - i ;
00520 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
00521         int ixsym = nm1 -  isym ;
00522 
00523         double fp = .5 * ( g[i] + g[isym] ) ;
00524         double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
00525 
00526         ff0[ix] = fp + fm ;
00527         ff0[ixsym] = fp - fm ;
00528             }
00529     
00530 //... cas particuliers:
00531         ff0[0] = g[0] - fmoins0 ;
00532         ff0[nm1] = g[0] + fmoins0 ;
00533         ff0[nm1s2] = g[nm1s2] ;
00534 
00535     }   // fin de la boucle sur theta 
00536 
00537 
00538 //--------------------------------------------------------------------
00539 //  partie sin(m phi) avec m impair : developpement en T_{2i}(x)
00540 //--------------------------------------------------------------------
00541 
00542     j++ ;
00543 
00544     if ( j != n1f-1 ) {  
00545 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00546 //  pas nuls 
00547 
00548     for (k=0; k<n2c; k++) {
00549 
00550         int i0 = n2n3c * j + n3c * k ; // indice de depart 
00551         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00552 
00553         i0 = n2n3f * j + n3f * k ; // indice de depart 
00554         double* ff0 = ff + i0 ;    // tableau resultat
00555 
00556 /*
00557  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00558  *     reliee a x par  x = cos(psi/2)   et F(psi) = f(x(psi)).  
00559  */
00560 
00561 // Calcul des coefficients de Fourier de la fonction 
00562 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00563 // en fonction des coefficients de Tchebyshev de f:
00564 
00565 // Coefficients impairs de G
00566 //--------------------------
00567  
00568         double c1 = cf0[1] ;
00569 
00570             double som = 0;
00571         ff0[1] = 0 ;
00572             for ( i = 3; i < nr; i += 2 ) {
00573             ff0[i] = cf0[i] - c1 ;
00574         som += ff0[i] ;
00575             }   
00576 
00577 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00578         double fmoins0 = nm1s2 * c1 + som ;
00579 
00580 // Coef. impairs de G
00581 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
00582 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00583             g[1] = 0 ;
00584             for ( i = 3; i < nr; i += 2 ) {
00585         g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
00586             }
00587             g[nr] = 0 ; 
00588 
00589 
00590 // Coefficients pairs de G
00591 //------------------------
00592 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00593 //   f.
00594 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
00595 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00596 
00597         g[0] = cf0[0] ;
00598             for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * cf0[i] ; 
00599             g[nm1] = cf0[nm1] ;
00600 
00601 // Transformation de Fourier inverse de G 
00602 //---------------------------------------
00603 
00604 // FFT inverse
00605             F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
00606 
00607 // Valeurs de f deduites de celles de G
00608 //-------------------------------------
00609 
00610             for ( i = 1; i < nm1s2 ; i++ ) {
00611 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
00612         int isym = nm1 - i ; 
00613 // ... indice (dans le tableau ff0) du point x correspondant a psi
00614         int ix = nm1 - i ;
00615 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
00616         int ixsym = nm1 -  isym ;
00617 
00618         double fp = .5 * ( g[i] + g[isym] ) ;
00619         double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
00620 
00621         ff0[ix] = fp + fm ;
00622         ff0[ixsym] = fp - fm ;
00623             }
00624     
00625 //... cas particuliers:
00626         ff0[0] = g[0] - fmoins0 ;
00627         ff0[nm1] = g[0] + fmoins0 ;
00628         ff0[nm1s2] = g[nm1s2] ;
00629 
00630     }   // fin de la boucle sur theta 
00631 
00632     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00633          // coef en phi n'etaient pas nuls)
00634 
00635 // On passe au cas m impair suivant:
00636     j+=3 ;
00637 
00638    }    // fin de la boucle sur les cas m impair
00639 
00640     // Menage
00641     free (t1) ;
00642     free (g) ;
00643     
00644 }
00645 

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