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 char cftcossinc_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinc.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $" ;
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.3 2014/10/13 08:53:19 j_novak Exp $
90  * $Log: cftcossinc.C,v $
91  * Revision 1.3 2014/10/13 08:53:19 j_novak
92  * Lorene classes and functions now belong to the namespace Lorene.
93  *
94  * Revision 1.2 2014/10/06 15:18:48 j_novak
95  * Modified #include directives to use c++ syntax.
96  *
97  * Revision 1.1 2004/12/21 17:06:02 j_novak
98  * Added all files for using fftw3.
99  *
100  * Revision 1.1 2004/11/23 15:13:50 m_forot
101  * Added the bases for the cases without any equatorial symmetry
102  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
103  *
104  *
105  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinc.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $
106  *
107  */
108 
109 // headers du C
110 #include <cstdlib>
111 #include <fftw3.h>
112 
113 //Lorene prototypes
114 #include "tbl.h"
115 
116 // Prototypage des sous-routines utilisees:
117 namespace Lorene {
118 fftw_plan prepare_fft(int, Tbl*&) ;
119 double* cheb_ini(const int) ;
120 //*****************************************************************************
121 
122 void cftcossinc(const int* deg, const int* dimf, double* ff, const int* dimc,
123  double* cf)
124 {
125 
126 int i, j, k ;
127 
128 // Dimensions des tableaux ff et cf :
129  int n1f = dimf[0] ;
130  int n2f = dimf[1] ;
131  int n3f = dimf[2] ;
132  int n1c = dimc[0] ;
133  int n2c = dimc[1] ;
134  int n3c = dimc[2] ;
135 
136 // Nombre de degres de liberte en theta :
137  int nt = deg[1] ;
138 
139 // Tests de dimension:
140  if (nt > n2f) {
141  cout << "cftcossinc: nt > n2f : nt = " << nt << " , n2f = "
142  << n2f << endl ;
143  abort () ;
144  exit(-1) ;
145  }
146  if (nt > n2c) {
147  cout << "cftcossinc: nt > n2c : nt = " << nt << " , n2c = "
148  << n2c << endl ;
149  abort () ;
150  exit(-1) ;
151  }
152  if (n1f > n1c) {
153  cout << "cftcossinc: n1f > n1c : n1f = " << n1f << " , n1c = "
154  << n1c << endl ;
155  abort () ;
156  exit(-1) ;
157  }
158  if (n3f > n3c) {
159  cout << "cftcossinc: n3f > n3c : n3f = " << n3f << " , n3c = "
160  << n3c << endl ;
161  abort () ;
162  exit(-1) ;
163  }
164 
165 // Nombre de points pour la FFT:
166  int nm1 = nt - 1;
167  int nm1s2 = nm1 / 2;
168 
169 // Recherche des tables pour la FFT:
170  Tbl* pg = 0x0 ;
171  fftw_plan p = prepare_fft(nm1, pg) ;
172  Tbl& g = *pg ;
173 
174 // Recherche de la table des sin(psi) :
175  double* sinp = cheb_ini(nt);
176 
177 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
178 // et 0 a dimf[2])
179 
180  int n2n3f = n2f * n3f ;
181  int n2n3c = n2c * n3c ;
182 
183 //=======================================================================
184 // Cas m pair
185 //=======================================================================
186 
187  j = 0 ;
188 
189  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
190  // (car nul)
191 
192 //------------------------------------------------------------------------
193 // partie cos(m phi) avec m pair : transformation en cos(l theta)
194 //------------------------------------------------------------------------
195 
196 
197  for (k=0; k<n3f; k++) {
198 
199  int i0 = n2n3f * j + k ; // indice de depart
200  double* ff0 = ff + i0 ; // tableau des donnees a transformer
201 
202  i0 = n2n3c * j + k ; // indice de depart
203  double* cf0 = cf + i0 ; // tableau resultat
204 
205 
206 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
207  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
208 
209 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
210 //---------------------------------------------
211  for ( i = 1; i < nm1s2 ; i++ ) {
212  int isym = nm1 - i ;
213  int ix = n3f * i ;
214  int ixsym = n3f * isym ;
215 // ... F+(theta)
216  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
217 // ... F_(theta) sin(psi)
218  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
219  g.set(i) = fp + fms ;
220  g.set(isym) = fp - fms ;
221  }
222 //... cas particuliers:
223  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
224  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
225 
226 // Developpement de G en series de Fourier par une FFT
227 //----------------------------------------------------
228 
229  fftw_execute(p) ;
230 
231 // Coefficients pairs du developmt. cos(l theta) de f
232 //----------------------------------------------------
233 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
234 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
235 // de fftw) :
236 
237  double fac = 2./double(nm1) ;
238  cf0[0] = g(0) / double(nm1) ;
239  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
240  cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
241 
242 // Coefficients impairs du developmt. en cos(l theta) de f
243 //---------------------------------------------------------
244 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
245 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
246 // (si fftw donnait reellement les coef. en sinus, il faudrait le
247 // remplacer par un -2.)
248  fac *= 2. ;
249  cf0[n3c] = 0 ;
250  double som = 0;
251  for ( i = 3; i < nt; i += 2 ) {
252  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
253  som += cf0[n3c*i] ;
254  }
255 
256 // 2. Calcul de c_1 :
257  double c1 = ( fmoins0 - som ) / nm1s2 ;
258 
259 // 3. Coef. c_k avec k impair:
260  cf0[n3c] = c1 ;
261  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
262 
263 
264  } // fin de la boucle sur r
265 
266 //--------------------------------------------------------------------
267 // partie sin(m phi) avec m pair : transformation en cos(l theta)
268 //--------------------------------------------------------------------
269 
270  j++ ;
271 
272  if ( (j != 1) && (j != n1f-1 ) ) {
273 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
274 // pas nuls
275 
276  for (k=0; k<n3f; k++) {
277 
278  int i0 = n2n3f * j + k ; // indice de depart
279  double* ff0 = ff + i0 ; // tableau des donnees a transformer
280 
281  i0 = n2n3c * j + k ; // indice de depart
282  double* cf0 = cf + i0 ; // tableau resultat
283 
284 
285 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
286  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
287 
288 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
289 //---------------------------------------------
290  for ( i = 1; i < nm1s2 ; i++ ) {
291  int isym = nm1 - i ;
292  int ix = n3f * i ;
293  int ixsym = n3f * isym ;
294 // ... F+(theta)
295  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
296 // ... F_(theta) sin(psi)
297  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
298  g.set(i) = fp + fms ;
299  g.set(isym) = fp - fms ;
300  }
301 //... cas particuliers:
302  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
303  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
304 
305 // Developpement de G en series de Fourier par une FFT
306 //----------------------------------------------------
307 
308  fftw_execute(p) ;
309 
310 // Coefficients pairs du developmt. cos(l theta) de f
311 //----------------------------------------------------
312 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
313 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
314 // de fftw) :
315 
316  double fac = 2./double(nm1) ;
317  cf0[0] = g(0) / double(nm1) ;
318  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
319  cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
320 
321 // Coefficients impairs du developmt. en cos(l theta) de f
322 //---------------------------------------------------------
323 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
324 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
325 // (si fftw donnait reellement les coef. en sinus, il faudrait le
326 // remplacer par un -2.)
327  fac *= 2. ;
328  cf0[n3c] = 0 ;
329  double som = 0;
330  for ( i = 3; i < nt; i += 2 ) {
331  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
332  som += cf0[n3c*i] ;
333  }
334 
335 // 2. Calcul de c_1 :
336  double c1 = ( fmoins0 - som ) / nm1s2 ;
337 
338 // 3. Coef. c_k avec k impair:
339  cf0[n3c] = c1 ;
340  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
341 
342 
343  } // fin de la boucle sur r
344  } // fin du cas ou le calcul etait necessaire (i.e. ou les
345  // coef en phi n'etaient pas nuls)
346 
347 // On passe au cas m pair suivant:
348  j+=3 ;
349 
350  } // fin de la boucle sur les cas m pair
351 
352  if (n1f<=3) // cas m=0 seulement (symetrie axiale)
353  return ;
354 
355 //=======================================================================
356 // Cas m impair
357 //=======================================================================
358 
359  j = 2 ;
360 
361  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
362  // (car nul)
363 
364 //--------------------------------------------------------------------
365 // partie cos(m phi) avec m impair : transformation en sin(l) theta)
366 //--------------------------------------------------------------------
367 
368  for (k=0; k<n3f; k++) {
369 
370  int i0 = n2n3f * j + k ; // indice de depart
371  double* ff0 = ff + i0 ; // tableau des donnees a transformer
372 
373  i0 = n2n3c * j + k ; // indice de depart
374  double* cf0 = cf + i0 ; // tableau resultat
375 
376 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
377 //---------------------------------------------
378 
379  for ( i = 1; i < nm1s2 ; i++ ) {
380  int isym = nm1 - i ;
381  int ix = n3f * i ;
382  int ixsym = n3f * isym ;
383 // ... F+(theta)
384  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
385 // ... F_(theta) sin(theta)
386  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
387  g.set(i) = fp + fms ;
388  g.set(isym) = fp - fms ;
389  }
390 //... cas particuliers:
391  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
392  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
393 
394 // Developpement de G en series de Fourier par une FFT
395 //----------------------------------------------------
396 
397  fftw_execute(p) ;
398 
399 // Coefficients pairs du developmt. sin(l theta) de f
400 //----------------------------------------------------
401 // Ces coefficients sont egaux aux coefficients en sinus du developpement
402 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
403 // de fftw) :
404 
405  double fac = -2. / double(nm1) ;
406  cf0[0] = 0. ;
407  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
408  cf0[n3c*nm1] = 0. ;
409 
410 // Coefficients impairs du developmt. en sin(l theta) de f
411 //---------------------------------------------------------
412 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
413 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
414 // (si fftw donnait reellement les coef. en sinus, il faudrait le
415 // remplacer par un -2.)
416 
417  cf0[n3c] = -fac * g(0);
418  fac *= -2. ;
419  for ( i = 3; i < nt; i += 2 ) {
420  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
421  }
422 
423  } // fin de la boucle sur r
424 
425 //------------------------------------------------------------------------
426 // partie sin(m phi) avec m impair : transformation en sin(l theta)
427 //------------------------------------------------------------------------
428 
429  j++ ;
430 
431  if ( j != n1f-1 ) {
432 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
433 // pas nuls
434 
435  for (k=0; k<n3f; k++) {
436 
437  int i0 = n2n3f * j + k ; // indice de depart
438  double* ff0 = ff + i0 ; // tableau des donnees a transformer
439 
440  i0 = n2n3c * j + k ; // indice de depart
441  double* cf0 = cf + i0 ; // tableau resultat
442 
443 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
444 //---------------------------------------------
445  for ( i = 1; i < nm1s2 ; i++ ) {
446  int isym = nm1 - i ;
447  int ix = n3f * i ;
448  int ixsym = n3f * isym ;
449 // ... F+(theta)
450  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
451 // ... F_(theta) sin(theta)
452  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
453  g.set(i) = fp + fms ;
454  g.set(isym) = fp - fms ;
455  }
456 //... cas particuliers:
457  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
458  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
459 
460 // Developpement de G en series de Fourier par une FFT
461 //----------------------------------------------------
462 
463  fftw_execute(p) ;
464 
465 // Coefficients pairs du developmt. sin(l theta) de f
466 //----------------------------------------------------
467 // Ces coefficients sont egaux aux coefficients en sinus du developpement
468 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
469 // de fftw) :
470 
471  double fac = -2. / double(nm1) ;
472  cf0[0] = 0. ;
473  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
474  cf0[n3c*nm1] = 0. ;
475 
476 // Coefficients impairs du developmt. en sin(l theta) de f
477 //---------------------------------------------------------
478 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
479 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
480 // (si fftw donnait reellement les coef. en sinus, il faudrait le
481 // remplacer par un -2.)
482 
483  cf0[n3c] = -fac * g(0);
484  fac *= -2. ;
485  for ( i = 3; i < nt; i += 2 ) {
486  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
487  }
488 
489  } // fin de la boucle sur r
490 
491  } // fin du cas ou le calcul etait necessaire (i.e. ou les
492  // coef en phi n'etaient pas nuls)
493 
494 
495 // On passe au cas m impair suivant:
496  j+=3 ;
497 
498  } // fin de la boucle sur les cas m impair
499 
500 }
501 }
Lorene
Lorene prototypes.
Definition: app_hor.h:64