LORENE
FFTW3/cftcossinc.C
1 /*
2  * Copyright (c) 1999-2002 Eric Gourgoulhon
3  *
4  * This file is part of LORENE.
5  *
6  * LORENE is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * LORENE is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with LORENE; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  *
20  */
21 
22 
23 
24 
25 /*
26  * Transformation en sin(l*theta) ou cos(l*theta) (suivant la parite
27  * de l'indice m en phi) sur le deuxieme indice (theta)
28  * d'un tableau 3-D representant une fonction sans symetrie par rapport
29  * au plan z=0.
30  * Utilise la bibliotheque fftw
31  *
32  * Entree:
33  * -------
34  * int* deg : tableau du nombre effectif de degres de liberte dans chacune
35  * des 3 dimensions: le nombre de points de collocation
36  * en theta est nt = deg[1] et doit etre de la forme
37  * nt = 2*p + 1
38  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
39  * dimensions.
40  * On doit avoir dimf[1] >= deg[1] = nt.
41  *
42  * double* ff : tableau des valeurs de la fonction aux nt points de
43  * de collocation
44  *
45  * theta_l = pi l/(nt-1) 0 <= l <= nt-1
46  *
47  * L'espace memoire correspondant a ce
48  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
49  * etre alloue avant l'appel a la routine.
50  * Les valeurs de la fonction doivent etre stokees
51  * dans le tableau ff comme suit
52  * f( theta_l ) = ff[ dimf[1]*dimf[2] * m + k + dimf[2] * l ]
53  * ou m et k sont les indices correspondant a
54  * phi et r respectivement.
55  * NB: cette routine suppose que la transformation en phi a deja ete
56  * effectuee: ainsi m est un indice de Fourier, non un indice de
57  * point de collocation en phi.
58  *
59  * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
60  * dimensions.
61  * On doit avoir dimc[1] >= deg[1] = nt.
62  * Sortie:
63  * -------
64  * double* cf : tableau des coefficients c_l de la fonction definis
65  * comme suit (a r et phi fixes)
66  *
67  * pour m pair:
68  * f(theta) = som_{l=0}^{nt-1} c_l cos(l theta ) .
69  * pour m impair:
70  * f(theta) = som_{l=0}^{nt-1} c_l sin( l theta ) .
71  *
72  * L'espace memoire correspondant a ce
73  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
74  * etre alloue avant l'appel a la routine.
75  * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
76  * le tableau cf comme suit
77  * c_l = cf[ dimc[1]*dimc[2] * m + k + dimc[2] * l ]
78  * ou m et k sont les indices correspondant a
79  * phi et r respectivement.
80  * Pour m pair, c_0 = c_{nt-1} = 0.
81  * Pour m impair, c_{nt-1} = 0.
82  *
83  * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
84  * seul tableau, qui constitue une entree/sortie.
85  *
86  */
87 
88 /*
89  * $Id: cftcossinc.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
90  * $Log: cftcossinc.C,v $
91  * Revision 1.4 2016/12/05 16:18:05 j_novak
92  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
93  *
94  * Revision 1.3 2014/10/13 08:53:19 j_novak
95  * Lorene classes and functions now belong to the namespace Lorene.
96  *
97  * Revision 1.2 2014/10/06 15:18:48 j_novak
98  * Modified #include directives to use c++ syntax.
99  *
100  * Revision 1.1 2004/12/21 17:06:02 j_novak
101  * Added all files for using fftw3.
102  *
103  * Revision 1.1 2004/11/23 15:13:50 m_forot
104  * Added the bases for the cases without any equatorial symmetry
105  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
106  *
107  *
108  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinc.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
109  *
110  */
111 
112 // headers du C
113 #include <cstdlib>
114 #include <fftw3.h>
115 
116 //Lorene prototypes
117 #include "tbl.h"
118 
119 // Prototypage des sous-routines utilisees:
120 namespace Lorene {
121 fftw_plan prepare_fft(int, Tbl*&) ;
122 double* cheb_ini(const int) ;
123 //*****************************************************************************
124 
125 void cftcossinc(const int* deg, const int* dimf, double* ff, const int* dimc,
126  double* cf)
127 {
128 
129 int i, j, k ;
130 
131 // Dimensions des tableaux ff et cf :
132  int n1f = dimf[0] ;
133  int n2f = dimf[1] ;
134  int n3f = dimf[2] ;
135  int n1c = dimc[0] ;
136  int n2c = dimc[1] ;
137  int n3c = dimc[2] ;
138 
139 // Nombre de degres de liberte en theta :
140  int nt = deg[1] ;
141 
142 // Tests de dimension:
143  if (nt > n2f) {
144  cout << "cftcossinc: nt > n2f : nt = " << nt << " , n2f = "
145  << n2f << endl ;
146  abort () ;
147  exit(-1) ;
148  }
149  if (nt > n2c) {
150  cout << "cftcossinc: nt > n2c : nt = " << nt << " , n2c = "
151  << n2c << endl ;
152  abort () ;
153  exit(-1) ;
154  }
155  if (n1f > n1c) {
156  cout << "cftcossinc: n1f > n1c : n1f = " << n1f << " , n1c = "
157  << n1c << endl ;
158  abort () ;
159  exit(-1) ;
160  }
161  if (n3f > n3c) {
162  cout << "cftcossinc: n3f > n3c : n3f = " << n3f << " , n3c = "
163  << n3c << endl ;
164  abort () ;
165  exit(-1) ;
166  }
167 
168 // Nombre de points pour la FFT:
169  int nm1 = nt - 1;
170  int nm1s2 = nm1 / 2;
171 
172 // Recherche des tables pour la FFT:
173  Tbl* pg = 0x0 ;
174  fftw_plan p = prepare_fft(nm1, pg) ;
175  Tbl& g = *pg ;
176 
177 // Recherche de la table des sin(psi) :
178  double* sinp = cheb_ini(nt);
179 
180 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
181 // et 0 a dimf[2])
182 
183  int n2n3f = n2f * n3f ;
184  int n2n3c = n2c * n3c ;
185 
186 //=======================================================================
187 // Cas m pair
188 //=======================================================================
189 
190  j = 0 ;
191 
192  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
193  // (car nul)
194 
195 //------------------------------------------------------------------------
196 // partie cos(m phi) avec m pair : transformation en cos(l theta)
197 //------------------------------------------------------------------------
198 
199 
200  for (k=0; k<n3f; k++) {
201 
202  int i0 = n2n3f * j + k ; // indice de depart
203  double* ff0 = ff + i0 ; // tableau des donnees a transformer
204 
205  i0 = n2n3c * j + k ; // indice de depart
206  double* cf0 = cf + i0 ; // tableau resultat
207 
208 
209 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
210  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
211 
212 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
213 //---------------------------------------------
214  for ( i = 1; i < nm1s2 ; i++ ) {
215  int isym = nm1 - i ;
216  int ix = n3f * i ;
217  int ixsym = n3f * isym ;
218 // ... F+(theta)
219  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
220 // ... F_(theta) sin(psi)
221  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
222  g.set(i) = fp + fms ;
223  g.set(isym) = fp - fms ;
224  }
225 //... cas particuliers:
226  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
227  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
228 
229 // Developpement de G en series de Fourier par une FFT
230 //----------------------------------------------------
231 
232  fftw_execute(p) ;
233 
234 // Coefficients pairs du developmt. cos(l theta) de f
235 //----------------------------------------------------
236 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
237 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
238 // de fftw) :
239 
240  double fac = 2./double(nm1) ;
241  cf0[0] = g(0) / double(nm1) ;
242  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
243  cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
244 
245 // Coefficients impairs du developmt. en cos(l theta) de f
246 //---------------------------------------------------------
247 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
248 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
249 // (si fftw donnait reellement les coef. en sinus, il faudrait le
250 // remplacer par un -2.)
251  fac *= 2. ;
252  cf0[n3c] = 0 ;
253  double som = 0;
254  for ( i = 3; i < nt; i += 2 ) {
255  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
256  som += cf0[n3c*i] ;
257  }
258 
259 // 2. Calcul de c_1 :
260  double c1 = ( fmoins0 - som ) / nm1s2 ;
261 
262 // 3. Coef. c_k avec k impair:
263  cf0[n3c] = c1 ;
264  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
265 
266 
267  } // fin de la boucle sur r
268 
269 //--------------------------------------------------------------------
270 // partie sin(m phi) avec m pair : transformation en cos(l theta)
271 //--------------------------------------------------------------------
272 
273  j++ ;
274 
275  if ( (j != 1) && (j != n1f-1 ) ) {
276 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
277 // pas nuls
278 
279  for (k=0; k<n3f; k++) {
280 
281  int i0 = n2n3f * j + k ; // indice de depart
282  double* ff0 = ff + i0 ; // tableau des donnees a transformer
283 
284  i0 = n2n3c * j + k ; // indice de depart
285  double* cf0 = cf + i0 ; // tableau resultat
286 
287 
288 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
289  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
290 
291 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
292 //---------------------------------------------
293  for ( i = 1; i < nm1s2 ; i++ ) {
294  int isym = nm1 - i ;
295  int ix = n3f * i ;
296  int ixsym = n3f * isym ;
297 // ... F+(theta)
298  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
299 // ... F_(theta) sin(psi)
300  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
301  g.set(i) = fp + fms ;
302  g.set(isym) = fp - fms ;
303  }
304 //... cas particuliers:
305  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
306  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
307 
308 // Developpement de G en series de Fourier par une FFT
309 //----------------------------------------------------
310 
311  fftw_execute(p) ;
312 
313 // Coefficients pairs du developmt. cos(l theta) de f
314 //----------------------------------------------------
315 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
316 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
317 // de fftw) :
318 
319  double fac = 2./double(nm1) ;
320  cf0[0] = g(0) / double(nm1) ;
321  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
322  cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
323 
324 // Coefficients impairs du developmt. en cos(l theta) de f
325 //---------------------------------------------------------
326 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
327 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
328 // (si fftw donnait reellement les coef. en sinus, il faudrait le
329 // remplacer par un -2.)
330  fac *= 2. ;
331  cf0[n3c] = 0 ;
332  double som = 0;
333  for ( i = 3; i < nt; i += 2 ) {
334  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
335  som += cf0[n3c*i] ;
336  }
337 
338 // 2. Calcul de c_1 :
339  double c1 = ( fmoins0 - som ) / nm1s2 ;
340 
341 // 3. Coef. c_k avec k impair:
342  cf0[n3c] = c1 ;
343  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
344 
345 
346  } // fin de la boucle sur r
347  } // fin du cas ou le calcul etait necessaire (i.e. ou les
348  // coef en phi n'etaient pas nuls)
349 
350 // On passe au cas m pair suivant:
351  j+=3 ;
352 
353  } // fin de la boucle sur les cas m pair
354 
355  if (n1f<=3) // cas m=0 seulement (symetrie axiale)
356  return ;
357 
358 //=======================================================================
359 // Cas m impair
360 //=======================================================================
361 
362  j = 2 ;
363 
364  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
365  // (car nul)
366 
367 //--------------------------------------------------------------------
368 // partie cos(m phi) avec m impair : transformation en sin(l) theta)
369 //--------------------------------------------------------------------
370 
371  for (k=0; k<n3f; k++) {
372 
373  int i0 = n2n3f * j + k ; // indice de depart
374  double* ff0 = ff + i0 ; // tableau des donnees a transformer
375 
376  i0 = n2n3c * j + k ; // indice de depart
377  double* cf0 = cf + i0 ; // tableau resultat
378 
379 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
380 //---------------------------------------------
381 
382  for ( i = 1; i < nm1s2 ; i++ ) {
383  int isym = nm1 - i ;
384  int ix = n3f * i ;
385  int ixsym = n3f * isym ;
386 // ... F+(theta)
387  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
388 // ... F_(theta) sin(theta)
389  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
390  g.set(i) = fp + fms ;
391  g.set(isym) = fp - fms ;
392  }
393 //... cas particuliers:
394  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
395  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
396 
397 // Developpement de G en series de Fourier par une FFT
398 //----------------------------------------------------
399 
400  fftw_execute(p) ;
401 
402 // Coefficients pairs du developmt. sin(l theta) de f
403 //----------------------------------------------------
404 // Ces coefficients sont egaux aux coefficients en sinus du developpement
405 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
406 // de fftw) :
407 
408  double fac = -2. / double(nm1) ;
409  cf0[0] = 0. ;
410  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
411  cf0[n3c*nm1] = 0. ;
412 
413 // Coefficients impairs du developmt. en sin(l theta) de f
414 //---------------------------------------------------------
415 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
416 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
417 // (si fftw donnait reellement les coef. en sinus, il faudrait le
418 // remplacer par un -2.)
419 
420  cf0[n3c] = -fac * g(0);
421  fac *= -2. ;
422  for ( i = 3; i < nt; i += 2 ) {
423  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
424  }
425 
426  } // fin de la boucle sur r
427 
428 //------------------------------------------------------------------------
429 // partie sin(m phi) avec m impair : transformation en sin(l theta)
430 //------------------------------------------------------------------------
431 
432  j++ ;
433 
434  if ( j != n1f-1 ) {
435 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
436 // pas nuls
437 
438  for (k=0; k<n3f; k++) {
439 
440  int i0 = n2n3f * j + k ; // indice de depart
441  double* ff0 = ff + i0 ; // tableau des donnees a transformer
442 
443  i0 = n2n3c * j + k ; // indice de depart
444  double* cf0 = cf + i0 ; // tableau resultat
445 
446 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
447 //---------------------------------------------
448  for ( i = 1; i < nm1s2 ; i++ ) {
449  int isym = nm1 - i ;
450  int ix = n3f * i ;
451  int ixsym = n3f * isym ;
452 // ... F+(theta)
453  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
454 // ... F_(theta) sin(theta)
455  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
456  g.set(i) = fp + fms ;
457  g.set(isym) = fp - fms ;
458  }
459 //... cas particuliers:
460  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
461  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
462 
463 // Developpement de G en series de Fourier par une FFT
464 //----------------------------------------------------
465 
466  fftw_execute(p) ;
467 
468 // Coefficients pairs du developmt. sin(l theta) de f
469 //----------------------------------------------------
470 // Ces coefficients sont egaux aux coefficients en sinus du developpement
471 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
472 // de fftw) :
473 
474  double fac = -2. / double(nm1) ;
475  cf0[0] = 0. ;
476  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
477  cf0[n3c*nm1] = 0. ;
478 
479 // Coefficients impairs du developmt. en sin(l theta) de f
480 //---------------------------------------------------------
481 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
482 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
483 // (si fftw donnait reellement les coef. en sinus, il faudrait le
484 // remplacer par un -2.)
485 
486  cf0[n3c] = -fac * g(0);
487  fac *= -2. ;
488  for ( i = 3; i < nt; i += 2 ) {
489  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
490  }
491 
492  } // fin de la boucle sur r
493 
494  } // fin du cas ou le calcul etait necessaire (i.e. ou les
495  // coef en phi n'etaient pas nuls)
496 
497 
498 // On passe au cas m impair suivant:
499  j+=3 ;
500 
501  } // fin de la boucle sur les cas m impair
502 
503 }
504 }
Lorene prototypes.
Definition: app_hor.h:67