citcossins.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 citcossins_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossins.C,v 1.1 2004/12/21 17:06:03 j_novak Exp $" ;
00024 
00025 
00026 /*
00027  * Transformation inverse cos(l*theta) ou sin(l*theta) (suivant la 
00028  *  parite 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 bibliotheuqe 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* dimc  : tableau du nombre d'elements de cc dans chacune des trois 
00040  *            dimensions.
00041  *        On doit avoir  dimc[1] >= deg[1] = nt. 
00042  *
00043  *   double* cf :   tableau des coefficients c_l de la fonction definis
00044  *            comme suit (a r et phi fixes)
00045  *
00046  *            pour m pair:
00047  *          f(theta) = som_{l=0}^{nt-1} c_l sin( l theta ) . 
00048  *            pour m impair:
00049  *          f(theta) = som_{l=0}^{nt-2} c_l cos( l theta ) . 
00050  *
00051  *            L'espace memoire correspondant a ce
00052  *                        pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit 
00053  *            avoir ete alloue avant l'appel a la routine.   
00054  *           Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
00055  *           le tableau cf comme suit
00056  *                c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
00057  *           ou j et k sont les indices correspondant a
00058  *           phi et r respectivement.
00059  *
00060  *   int* dimf  : tableau du nombre d'elements de ff dans chacune des trois 
00061  *            dimensions.
00062  *        On doit avoir  dimf[1] >= deg[1] = nt. 
00063  *
00064  * Sortie:
00065  * -------
00066  *   double* ff : tableau des valeurs de la fonction aux nt points de
00067  *                        de collocation
00068  *
00069  *            theta_l =  pi l/(nt-1)       0 <= l <= nt-1 
00070  *
00071  *            L'espace memoire correspondant a ce
00072  *                        pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit 
00073  *            avoir ete alloue avant l'appel a la routine.   
00074  *            Les valeurs de la fonction sont stokees
00075  *            dans le tableau ff comme suit
00076  *          f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
00077  *           ou j et k sont les indices correspondant a
00078  *           phi et r respectivement.
00079  *
00080  * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un 
00081  *     seul tableau, qui constitue une entree/sortie.
00082  *
00083  */
00084 
00085 /*
00086  * $Id: citcossins.C,v 1.1 2004/12/21 17:06:03 j_novak Exp $
00087  * $Log: citcossins.C,v $
00088  * Revision 1.1  2004/12/21 17:06:03  j_novak
00089  * Added all files for using fftw3.
00090  *
00091  * Revision 1.1  2004/11/23 15:13:50  m_forot
00092  * Added the bases for the cases without any equatorial symmetry
00093  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
00094  *
00095  *
00096  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossins.C,v 1.1 2004/12/21 17:06:03 j_novak Exp $
00097  *
00098  */
00099 // headers du C
00100 #include <stdlib.h>
00101 #include <fftw3.h>
00102 
00103 //Lorene prototypes
00104 #include "tbl.h"
00105 
00106 // Prototypage des sous-routines utilisees:
00107 fftw_plan back_fft(int, Tbl*&) ;
00108 double* cheb_ini(const int) ;
00109 //*****************************************************************************
00110 
00111 void citcossins(const int* deg, const int* dimc, double* cf, const int* dimf,
00112            double* ff)
00113 {
00114 
00115 int i, j, k ;
00116 
00117 // Dimensions des tableaux ff et cf  :
00118     int n1f = dimf[0] ;
00119     int n2f = dimf[1] ;
00120     int n3f = dimf[2] ;
00121     int n1c = dimc[0] ;
00122     int n2c = dimc[1] ;
00123     int n3c = dimc[2] ;
00124 
00125 // Nombres de degres de liberte en theta :    
00126     int nt = deg[1] ;
00127     
00128 // Tests de dimension:
00129     if (nt > n2f) {
00130     cout << "citcossins: nt > n2f : nt = " << nt << " ,  n2f = " 
00131     << n2f << endl ;
00132     abort () ;
00133     exit(-1) ;
00134     }
00135     if (nt > n2c) {
00136     cout << "citcossins: nt > n2c : nt = " << nt << " ,  n2c = " 
00137     << n2c << endl ;
00138     abort () ;
00139     exit(-1) ;
00140     }
00141     if (n1c > n1f) {
00142     cout << "citcossins: n1c > n1f : n1c = " << n1c << " ,  n1f = " 
00143     << n1f << endl ;
00144     abort () ;
00145     exit(-1) ;
00146     }
00147     if (n3c > n3f) {
00148     cout << "citcossins: n3c > n3f : n3c = " << n3c << " ,  n3f = " 
00149     << n3f << endl ;
00150     abort () ;
00151     exit(-1) ;
00152     }
00153 
00154 // Nombre de points pour la FFT:
00155     int nm1 = nt - 1;
00156     int nm1s2 = nm1 / 2;
00157 
00158 // Recherche des tables pour la FFT:
00159     Tbl* pg = 0x0 ;
00160     fftw_plan p = back_fft(nm1, pg) ;
00161     Tbl& g = *pg ;
00162 
00163 // Recherche de la table des sin(psi) :
00164     double* sinp = cheb_ini(nt);    
00165     
00166 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
00167 //           et 0 a dimf[2])
00168 
00169     int n2n3f = n2f * n3f ;
00170     int n2n3c = n2c * n3c ;
00171 
00172 //=======================================================================
00173 //              Cas m pair
00174 //=======================================================================
00175 
00176     j = 0 ;
00177     
00178     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00179             // (car nul)
00180 
00181 //--------------------------------------------------------------------------
00182 //  partie cos(m phi) avec m pair : transformation sin(l theta) inv.
00183 //--------------------------------------------------------------------------
00184 
00185         for (k=0; k<n3c; k++) {
00186 
00187         int i0 = n2n3c * j + k ; // indice de depart 
00188         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00189 
00190         i0 = n2n3f * j + k ; // indice de depart 
00191         double* ff0 = ff + i0 ;    // tableau resultat
00192          
00193 // Coefficients impairs de G
00194 //--------------------------
00195  
00196         for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
00197 
00198 // Coefficients pairs de G
00199 //------------------------
00200 
00201         g.set(0) = .5 * cf0[n3c] ;
00202             for ( i = 3; i < nt; i += 2 ) {
00203         g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
00204             }
00205         g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
00206         
00207 // Transformation de Fourier inverse de G 
00208 //---------------------------------------
00209 
00210 // FFT inverse
00211         fftw_execute(p) ;
00212 
00213 // Valeurs de f deduites de celles de G
00214 //-------------------------------------
00215 
00216             for ( i = 1; i < nm1s2 ; i++ ) {
00217 // ... indice du pt symetrique de psi par rapport a pi/2:
00218         int isym = nm1 - i ; 
00219     
00220         double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i]  ;
00221         double fm = 0.5 * ( g(i) - g(isym) ) ;
00222         ff0[ n3f*i ] = fp + fm ;
00223         ff0[ n3f*isym ] = fp - fm ;
00224             }
00225     
00226 //... cas particuliers:
00227         ff0[0] = 0. ;
00228         ff0[ n3f*nm1 ] = -2*g(0) ;
00229         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00230     }   // fin de la boucle sur r 
00231 
00232 //--------------------------------------------------------------------------
00233 //  partie sin(m phi) avec m pair : transformation sin(l theta) inv.
00234 //--------------------------------------------------------------------------
00235 
00236     j++ ;
00237 
00238     if ( j != n1f-1  ) {  
00239 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00240 //  pas nuls 
00241       
00242         for (k=0; k<n3c; k++) {
00243 
00244         int i0 = n2n3c * j + k ; // indice de depart 
00245         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00246 
00247         i0 = n2n3f * j + k ; // indice de depart 
00248         double* ff0 = ff + i0 ;    // tableau resultat
00249          
00250 // Coefficients impairs de G
00251 //--------------------------
00252  
00253         for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
00254 
00255 // Coefficients pairs de G
00256 //------------------------
00257 
00258         g.set(0) = .5 * cf0[n3c] ;
00259             for ( i = 3; i < nt; i += 2 ) {
00260         g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
00261             }
00262         g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
00263         
00264 // Transformation de Fourier inverse de G 
00265 //---------------------------------------
00266 
00267 // FFT inverse
00268         fftw_execute(p) ;
00269 
00270 // Valeurs de f deduites de celles de G
00271 //-------------------------------------
00272 
00273             for ( i = 1; i < nm1s2 ; i++ ) {
00274 // ... indice du pt symetrique de psi par rapport a pi/2:
00275         int isym = nm1 - i ; 
00276     
00277         double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i]  ;
00278         double fm = 0.5 * ( g(i) - g(isym) ) ;
00279         ff0[ n3f*i ] = fp + fm ;
00280         ff0[ n3f*isym ] = fp - fm ;
00281             }
00282     
00283 //... cas particuliers:
00284         ff0[0] = 0. ;
00285         ff0[ n3f*nm1 ] = -2*g(0) ;
00286         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00287     }   // fin de la boucle sur r 
00288 
00289     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00290          // coef en phi n'etaient pas nuls)
00291 
00292 // On passe au cas m pair suivant:
00293     j+=3 ;
00294 
00295    }    // fin de la boucle sur les cas m pair
00296 
00297     if (n1f<=3)         // cas m=0 seulement (symetrie axiale)
00298     return ;
00299     
00300 //=======================================================================
00301 //              Cas m impair
00302 //=======================================================================
00303 
00304     j = 2 ;
00305     
00306     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00307             // (car nul)
00308 
00309 //-----------------------------------------------------------------------
00310 //  partie cos(m phi) avec m impair : transformation cos( l theta) inverse
00311 //-----------------------------------------------------------------------
00312 
00313     for (k=0; k<n3c; k++) {
00314 
00315         int i0 = n2n3c * j + k ; // indice de depart 
00316         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00317 
00318         i0 = n2n3f * j + k ; // indice de depart 
00319         double* ff0 = ff + i0 ;    // tableau resultat
00320          
00321  
00322 // Coefficients impairs de G
00323 //--------------------------
00324  
00325         double c1 = cf0[n3c] ;
00326 
00327             double som = 0;
00328         ff0[n3f] = 0 ;
00329             for ( i = 3; i < nt; i += 2 ) {
00330             ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
00331         som += ff0[ n3f*i ] ;
00332             }   
00333 
00334 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
00335         double fmoins0 = nm1s2 * c1 + som ;
00336 
00337 // Coef. impairs de G
00338 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
00339 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00340             for ( i = 3; i < nt; i += 2 ) {
00341         g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
00342             }
00343 
00344 
00345 // Coefficients pairs de G
00346 //------------------------
00347 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00348 //   f.
00349 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
00350 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00351 
00352         g.set(0) = cf0[0] ;
00353             for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;  
00354             g.set(nm1s2) = cf0[ n3c*nm1 ] ;
00355 
00356 // Transformation de Fourier inverse de G 
00357 //---------------------------------------
00358 
00359 // FFT inverse
00360         fftw_execute(p) ;
00361 
00362 // Valeurs de f deduites de celles de G
00363 //-------------------------------------
00364 
00365             for ( i = 1; i < nm1s2 ; i++ ) {
00366 // ... indice du pt symetrique de psi par rapport a pi/2:
00367         int isym = nm1 - i ; 
00368     
00369         double fp = 0.5 * ( g(i) + g(isym) ) ;
00370         double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
00371         ff0[ n3f*i ] = fp + fm ;
00372         ff0[ n3f*isym ] = fp - fm ;
00373             }
00374     
00375 //... cas particuliers:
00376         ff0[0] = g(0) + fmoins0 ;
00377         ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
00378         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00379 
00380     }   // fin de la boucle sur r 
00381 
00382 //-----------------------------------------------------------------------
00383 //  partie sin(m phi) avec m impair : transformation cos(l theta) inverse
00384 //-----------------------------------------------------------------------
00385 
00386     j++ ;
00387 
00388     if ( (j != 1) && (j != n1f-1 ) ) {  
00389 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00390 //  pas nuls 
00391 
00392     for (k=0; k<n3c; k++) {
00393 
00394         int i0 = n2n3c * j + k ; // indice de depart 
00395         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00396 
00397         i0 = n2n3f * j + k ; // indice de depart 
00398         double* ff0 = ff + i0 ;    // tableau resultat
00399          
00400 // Coefficients impairs de G
00401 //--------------------------
00402  
00403         double c1 = cf0[n3c] ;
00404 
00405             double som = 0;
00406         ff0[n3f] = 0 ;
00407             for ( i = 3; i < nt; i += 2 ) {
00408             ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
00409         som += ff0[ n3f*i ] ;
00410             }   
00411 
00412 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
00413         double fmoins0 = nm1s2 * c1 + som ;
00414 
00415 // Coef. impairs de G
00416 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
00417 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00418             for ( i = 3; i < nt; i += 2 ) {
00419         g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
00420             }
00421 
00422 
00423 // Coefficients pairs de G
00424 //------------------------
00425 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00426 //   f.
00427 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
00428 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00429 
00430         g.set(0) = cf0[0] ;
00431             for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;  
00432             g.set(nm1s2) = cf0[ n3c*nm1 ] ;
00433 
00434 // Transformation de Fourier inverse de G 
00435 //---------------------------------------
00436 
00437 // FFT inverse
00438         fftw_execute(p) ;
00439 
00440 // Valeurs de f deduites de celles de G
00441 //-------------------------------------
00442 
00443             for ( i = 1; i < nm1s2 ; i++ ) {
00444 // ... indice du pt symetrique de psi par rapport a pi/2:
00445         int isym = nm1 - i ; 
00446     
00447         double fp = 0.5 * ( g(i) + g(isym) ) ;
00448         double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
00449         ff0[ n3f*i ] = fp + fm ;
00450         ff0[ n3f*isym ] = fp - fm ;
00451             }
00452     
00453 //... cas particuliers:
00454         ff0[0] = g(0) + fmoins0 ;
00455         ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
00456         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00457     }   // fin de la boucle sur r 
00458 
00459     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00460          // coef en phi n'etaient pas nuls)
00461 
00462 // On passe au cas m impair suivant:
00463     j+=3 ;
00464 
00465    }    // fin de la boucle sur les cas m impair
00466 
00467 }

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