citcossinc.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 citcossinc_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinc.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 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* 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 cos( l theta ) . 
00048  *            pour m impair:
00049  *          f(theta) = som_{l=0}^{nt-2} c_l sin( 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: citcossinc.C,v 1.1 2004/12/21 17:06:03 j_novak Exp $
00087  * $Log: citcossinc.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/citcossinc.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 citcossinc(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 << "citcossinc: nt > n2f : nt = " << nt << " ,  n2f = " 
00131     << n2f << endl ;
00132     abort () ;
00133     exit(-1) ;
00134     }
00135     if (nt > n2c) {
00136     cout << "citcossinc: nt > n2c : nt = " << nt << " ,  n2c = " 
00137     << n2c << endl ;
00138     abort () ;
00139     exit(-1) ;
00140     }
00141     if (n1c > n1f) {
00142     cout << "citcossinc: n1c > n1f : n1c = " << n1c << " ,  n1f = " 
00143     << n1f << endl ;
00144     abort () ;
00145     exit(-1) ;
00146     }
00147     if (n3c > n3f) {
00148     cout << "citcossinc: 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 cos( l theta) inverse
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  
00194 // Coefficients impairs de G
00195 //--------------------------
00196  
00197         double c1 = cf0[n3c] ;
00198 
00199             double som = 0;
00200         ff0[n3f] = 0 ;
00201             for ( i = 3; i < nt; i += 2 ) {
00202             ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
00203         som += ff0[ n3f*i ] ;
00204             }   
00205 
00206 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
00207         double fmoins0 = nm1s2 * c1 + som ;
00208 
00209 // Coef. impairs de G
00210 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
00211 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00212             for ( i = 3; i < nt; i += 2 ) {
00213         g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
00214             }
00215 
00216 
00217 // Coefficients pairs de G
00218 //------------------------
00219 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00220 //   f.
00221 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
00222 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00223 
00224         g.set(0) = cf0[0] ;
00225             for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;  
00226             g.set(nm1s2) = cf0[ n3c*nm1 ] ;
00227 
00228 // Transformation de Fourier inverse de G 
00229 //---------------------------------------
00230 
00231 // FFT inverse
00232         fftw_execute(p) ;
00233 
00234 // Valeurs de f deduites de celles de G
00235 //-------------------------------------
00236 
00237             for ( i = 1; i < nm1s2 ; i++ ) {
00238 // ... indice du pt symetrique de psi par rapport a pi/2:
00239         int isym = nm1 - i ; 
00240     
00241         double fp = 0.5 * ( g(i) + g(isym) ) ;
00242         double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
00243         ff0[ n3f*i ] = fp + fm ;
00244         ff0[ n3f*isym ] = fp - fm ;
00245             }
00246     
00247 //... cas particuliers:
00248         ff0[0] = g(0) + fmoins0 ;
00249         ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
00250         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00251     }   // fin de la boucle sur r 
00252 
00253 //-----------------------------------------------------------------------
00254 //  partie sin(m phi) avec m pair : transformation cos(l theta) inverse
00255 //-----------------------------------------------------------------------
00256 
00257     j++ ;
00258 
00259     if ( (j != 1) && (j != n1f-1 ) ) {  
00260 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00261 //  pas nuls 
00262 
00263     for (k=0; k<n3c; k++) {
00264 
00265         int i0 = n2n3c * j + k ; // indice de depart 
00266         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00267 
00268         i0 = n2n3f * j + k ; // indice de depart 
00269         double* ff0 = ff + i0 ;    // tableau resultat
00270          
00271 // Coefficients impairs de G
00272 //--------------------------
00273  
00274         double c1 = cf0[n3c] ;
00275 
00276             double som = 0;
00277         ff0[n3f] = 0 ;
00278             for ( i = 3; i < nt; i += 2 ) {
00279             ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
00280         som += ff0[ n3f*i ] ;
00281             }   
00282 
00283 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
00284         double fmoins0 = nm1s2 * c1 + som ;
00285 
00286 // Coef. impairs de G
00287 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
00288 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00289             for ( i = 3; i < nt; i += 2 ) {
00290         g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
00291             }
00292 
00293 
00294 // Coefficients pairs de G
00295 //------------------------
00296 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00297 //   f.
00298 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
00299 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00300 
00301         g.set(0) = cf0[0] ;
00302             for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;  
00303             g.set(nm1s2) = cf0[ n3c*nm1 ] ;
00304 
00305 // Transformation de Fourier inverse de G 
00306 //---------------------------------------
00307 
00308 // FFT inverse
00309         fftw_execute(p) ;
00310 
00311 // Valeurs de f deduites de celles de G
00312 //-------------------------------------
00313 
00314             for ( i = 1; i < nm1s2 ; i++ ) {
00315 // ... indice du pt symetrique de psi par rapport a pi/2:
00316         int isym = nm1 - i ; 
00317     
00318         double fp = 0.5 * ( g(i) + g(isym) ) ;
00319         double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
00320         ff0[ n3f*i ] = fp + fm ;
00321         ff0[ n3f*isym ] = fp - fm ;
00322             }
00323     
00324 //... cas particuliers:
00325         ff0[0] = g(0) + fmoins0 ;
00326         ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
00327         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00328     }   // fin de la boucle sur r 
00329 
00330     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00331          // coef en phi n'etaient pas nuls)
00332 
00333 // On passe au cas m pair suivant:
00334     j+=3 ;
00335 
00336    }    // fin de la boucle sur les cas m pair
00337 
00338 //##
00339     if (n1f<=3)         // cas m=0 seulement (symetrie axiale)
00340     return ;
00341     
00342 //=======================================================================
00343 //              Cas m impair
00344 //=======================================================================
00345 
00346     j = 2 ;
00347     
00348     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00349             // (car nul)
00350 
00351 //--------------------------------------------------------------------------
00352 //  partie cos(m phi) avec m impair : transformation sin(l theta) inv.
00353 //--------------------------------------------------------------------------
00354 
00355         for (k=0; k<n3c; k++) {
00356 
00357         int i0 = n2n3c * j + k ; // indice de depart 
00358         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00359 
00360         i0 = n2n3f * j + k ; // indice de depart 
00361         double* ff0 = ff + i0 ;    // tableau resultat
00362          
00363 // Coefficients impairs de G
00364 //--------------------------
00365  
00366         for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
00367 
00368 // Coefficients pairs de G
00369 //------------------------
00370 
00371         g.set(0) = .5 * cf0[n3c] ;
00372             for ( i = 3; i < nt; i += 2 ) {
00373         g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
00374             }
00375         g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
00376         
00377 // Transformation de Fourier inverse de G 
00378 //---------------------------------------
00379 
00380 // FFT inverse
00381         fftw_execute(p) ;
00382 
00383 // Valeurs de f deduites de celles de G
00384 //-------------------------------------
00385 
00386             for ( i = 1; i < nm1s2 ; i++ ) {
00387 // ... indice du pt symetrique de psi par rapport a pi/2:
00388         int isym = nm1 - i ; 
00389     
00390         double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i]  ;
00391         double fm = 0.5 * ( g(i) - g(isym) ) ;
00392         ff0[ n3f*i ] = fp + fm ;
00393         ff0[ n3f*isym ] = fp - fm ;
00394             }
00395     
00396 //... cas particuliers:
00397         ff0[0] = 0. ;
00398         ff0[ n3f*nm1 ] = -2*g(0) ;
00399         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00400     }   // fin de la boucle sur r 
00401 
00402 //--------------------------------------------------------------------------
00403 //  partie sin(m phi) avec m impair : transformation sin(l theta) inv.
00404 //--------------------------------------------------------------------------
00405 
00406     j++ ;
00407 
00408     if ( j != n1f-1  ) {  
00409 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00410 //  pas nuls 
00411       
00412         for (k=0; k<n3c; k++) {
00413 
00414         int i0 = n2n3c * j + k ; // indice de depart 
00415         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00416 
00417         i0 = n2n3f * j + k ; // indice de depart 
00418         double* ff0 = ff + i0 ;    // tableau resultat
00419          
00420 // Coefficients impairs de G
00421 //--------------------------
00422  
00423         for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
00424 
00425 // Coefficients pairs de G
00426 //------------------------
00427 
00428         g.set(0) = .5 * cf0[n3c] ;
00429             for ( i = 3; i < nt; i += 2 ) {
00430         g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
00431             }
00432         g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
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) ) / sinp[i]  ;
00448         double fm = 0.5 * ( g(i) - g(isym) ) ;
00449         ff0[ n3f*i ] = fp + fm ;
00450         ff0[ n3f*isym ] = fp - fm ;
00451             }
00452     
00453 //... cas particuliers:
00454         ff0[0] = 0. ;
00455         ff0[ n3f*nm1 ] = -2*g(0) ;
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