LORENE
FFTW3/citcossins.C
1 /*
2  * Copyright (c) 1999-2001 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 /*
27  * Transformation inverse cos(l*theta) ou sin(l*theta) (suivant la
28  * parite de l'indice m en phi) sur le deuxieme indice (theta)
29  * d'un tableau 3-D representant une fonction symetrique par rapport
30  * au plan z=0.
31  * Utilise la bibliotheuqe fftw.
32  *
33  * Entree:
34  * -------
35  * int* deg : tableau du nombre effectif de degres de liberte dans chacune
36  * des 3 dimensions: le nombre de points de collocation
37  * en theta est nt = deg[1] et doit etre de la forme
38  * nt = 2*p + 1
39  * int* dimc : tableau du nombre d'elements de cc dans chacune des trois
40  * dimensions.
41  * On doit avoir dimc[1] >= deg[1] = nt.
42  *
43  * double* cf : tableau des coefficients c_l de la fonction definis
44  * comme suit (a r et phi fixes)
45  *
46  * pour m pair:
47  * f(theta) = som_{l=0}^{nt-1} c_l sin( l theta ) .
48  * pour m impair:
49  * f(theta) = som_{l=0}^{nt-2} c_l cos( l theta ) .
50  *
51  * L'espace memoire correspondant a ce
52  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
53  * avoir ete alloue avant l'appel a la routine.
54  * Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
55  * le tableau cf comme suit
56  * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
57  * ou j et k sont les indices correspondant a
58  * phi et r respectivement.
59  *
60  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
61  * dimensions.
62  * On doit avoir dimf[1] >= deg[1] = nt.
63  *
64  * Sortie:
65  * -------
66  * double* ff : tableau des valeurs de la fonction aux nt points de
67  * de collocation
68  *
69  * theta_l = pi l/(nt-1) 0 <= l <= nt-1
70  *
71  * L'espace memoire correspondant a ce
72  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
73  * avoir ete alloue avant l'appel a la routine.
74  * Les valeurs de la fonction sont stokees
75  * dans le tableau ff comme suit
76  * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
77  * ou j et k sont les indices correspondant a
78  * phi et r respectivement.
79  *
80  * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
81  * seul tableau, qui constitue une entree/sortie.
82  *
83  */
84 
85 /*
86  * $Id: citcossins.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
87  * $Log: citcossins.C,v $
88  * Revision 1.4 2016/12/05 16:18:05 j_novak
89  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
90  *
91  * Revision 1.3 2014/10/13 08:53:20 j_novak
92  * Lorene classes and functions now belong to the namespace Lorene.
93  *
94  * Revision 1.2 2014/10/06 15:18:50 j_novak
95  * Modified #include directives to use c++ syntax.
96  *
97  * Revision 1.1 2004/12/21 17:06:03 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/citcossins.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
106  *
107  */
108 // headers du C
109 #include <cstdlib>
110 #include <fftw3.h>
111 
112 //Lorene prototypes
113 #include "tbl.h"
114 
115 // Prototypage des sous-routines utilisees:
116 namespace Lorene {
117 fftw_plan back_fft(int, Tbl*&) ;
118 double* cheb_ini(const int) ;
119 //*****************************************************************************
120 
121 void citcossins(const int* deg, const int* dimc, double* cf, const int* dimf,
122  double* ff)
123 {
124 
125 int i, j, k ;
126 
127 // Dimensions des tableaux ff et cf :
128  int n1f = dimf[0] ;
129  int n2f = dimf[1] ;
130  int n3f = dimf[2] ;
131  int n1c = dimc[0] ;
132  int n2c = dimc[1] ;
133  int n3c = dimc[2] ;
134 
135 // Nombres de degres de liberte en theta :
136  int nt = deg[1] ;
137 
138 // Tests de dimension:
139  if (nt > n2f) {
140  cout << "citcossins: nt > n2f : nt = " << nt << " , n2f = "
141  << n2f << endl ;
142  abort () ;
143  exit(-1) ;
144  }
145  if (nt > n2c) {
146  cout << "citcossins: nt > n2c : nt = " << nt << " , n2c = "
147  << n2c << endl ;
148  abort () ;
149  exit(-1) ;
150  }
151  if (n1c > n1f) {
152  cout << "citcossins: n1c > n1f : n1c = " << n1c << " , n1f = "
153  << n1f << endl ;
154  abort () ;
155  exit(-1) ;
156  }
157  if (n3c > n3f) {
158  cout << "citcossins: n3c > n3f : n3c = " << n3c << " , n3f = "
159  << n3f << endl ;
160  abort () ;
161  exit(-1) ;
162  }
163 
164 // Nombre de points pour la FFT:
165  int nm1 = nt - 1;
166  int nm1s2 = nm1 / 2;
167 
168 // Recherche des tables pour la FFT:
169  Tbl* pg = 0x0 ;
170  fftw_plan p = back_fft(nm1, pg) ;
171  Tbl& g = *pg ;
172 
173 // Recherche de la table des sin(psi) :
174  double* sinp = cheb_ini(nt);
175 
176 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
177 // et 0 a dimf[2])
178 
179  int n2n3f = n2f * n3f ;
180  int n2n3c = n2c * n3c ;
181 
182 //=======================================================================
183 // Cas m pair
184 //=======================================================================
185 
186  j = 0 ;
187 
188  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
189  // (car nul)
190 
191 //--------------------------------------------------------------------------
192 // partie cos(m phi) avec m pair : transformation sin(l theta) inv.
193 //--------------------------------------------------------------------------
194 
195  for (k=0; k<n3c; k++) {
196 
197  int i0 = n2n3c * j + k ; // indice de depart
198  double* cf0 = cf + i0 ; // tableau des donnees a transformer
199 
200  i0 = n2n3f * j + k ; // indice de depart
201  double* ff0 = ff + i0 ; // tableau resultat
202 
203 // Coefficients impairs de G
204 //--------------------------
205 
206  for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
207 
208 // Coefficients pairs de G
209 //------------------------
210 
211  g.set(0) = .5 * cf0[n3c] ;
212  for ( i = 3; i < nt; i += 2 ) {
213  g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
214  }
215  g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
216 
217 // Transformation de Fourier inverse de G
218 //---------------------------------------
219 
220 // FFT inverse
221  fftw_execute(p) ;
222 
223 // Valeurs de f deduites de celles de G
224 //-------------------------------------
225 
226  for ( i = 1; i < nm1s2 ; i++ ) {
227 // ... indice du pt symetrique de psi par rapport a pi/2:
228  int isym = nm1 - i ;
229 
230  double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
231  double fm = 0.5 * ( g(i) - g(isym) ) ;
232  ff0[ n3f*i ] = fp + fm ;
233  ff0[ n3f*isym ] = fp - fm ;
234  }
235 
236 //... cas particuliers:
237  ff0[0] = 0. ;
238  ff0[ n3f*nm1 ] = -2*g(0) ;
239  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
240  } // fin de la boucle sur r
241 
242 //--------------------------------------------------------------------------
243 // partie sin(m phi) avec m pair : transformation sin(l theta) inv.
244 //--------------------------------------------------------------------------
245 
246  j++ ;
247 
248  if ( j != n1f-1 ) {
249 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
250 // pas nuls
251 
252  for (k=0; k<n3c; k++) {
253 
254  int i0 = n2n3c * j + k ; // indice de depart
255  double* cf0 = cf + i0 ; // tableau des donnees a transformer
256 
257  i0 = n2n3f * j + k ; // indice de depart
258  double* ff0 = ff + i0 ; // tableau resultat
259 
260 // Coefficients impairs de G
261 //--------------------------
262 
263  for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
264 
265 // Coefficients pairs de G
266 //------------------------
267 
268  g.set(0) = .5 * cf0[n3c] ;
269  for ( i = 3; i < nt; i += 2 ) {
270  g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
271  }
272  g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
273 
274 // Transformation de Fourier inverse de G
275 //---------------------------------------
276 
277 // FFT inverse
278  fftw_execute(p) ;
279 
280 // Valeurs de f deduites de celles de G
281 //-------------------------------------
282 
283  for ( i = 1; i < nm1s2 ; i++ ) {
284 // ... indice du pt symetrique de psi par rapport a pi/2:
285  int isym = nm1 - i ;
286 
287  double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
288  double fm = 0.5 * ( g(i) - g(isym) ) ;
289  ff0[ n3f*i ] = fp + fm ;
290  ff0[ n3f*isym ] = fp - fm ;
291  }
292 
293 //... cas particuliers:
294  ff0[0] = 0. ;
295  ff0[ n3f*nm1 ] = -2*g(0) ;
296  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
297  } // fin de la boucle sur r
298 
299  } // fin du cas ou le calcul etait necessaire (i.e. ou les
300  // coef en phi n'etaient pas nuls)
301 
302 // On passe au cas m pair suivant:
303  j+=3 ;
304 
305  } // fin de la boucle sur les cas m pair
306 
307  if (n1f<=3) // cas m=0 seulement (symetrie axiale)
308  return ;
309 
310 //=======================================================================
311 // Cas m impair
312 //=======================================================================
313 
314  j = 2 ;
315 
316  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
317  // (car nul)
318 
319 //-----------------------------------------------------------------------
320 // partie cos(m phi) avec m impair : transformation cos( l theta) inverse
321 //-----------------------------------------------------------------------
322 
323  for (k=0; k<n3c; k++) {
324 
325  int i0 = n2n3c * j + k ; // indice de depart
326  double* cf0 = cf + i0 ; // tableau des donnees a transformer
327 
328  i0 = n2n3f * j + k ; // indice de depart
329  double* ff0 = ff + i0 ; // tableau resultat
330 
331 
332 // Coefficients impairs de G
333 //--------------------------
334 
335  double c1 = cf0[n3c] ;
336 
337  double som = 0;
338  ff0[n3f] = 0 ;
339  for ( i = 3; i < nt; i += 2 ) {
340  ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
341  som += ff0[ n3f*i ] ;
342  }
343 
344 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
345  double fmoins0 = nm1s2 * c1 + som ;
346 
347 // Coef. impairs de G
348 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
349 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
350  for ( i = 3; i < nt; i += 2 ) {
351  g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
352  }
353 
354 
355 // Coefficients pairs de G
356 //------------------------
357 // Ces coefficients sont egaux aux coefficients pairs du developpement de
358 // f.
359 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
360 // donnait exactement les coef. des cosinus, ce facteur serait 1.
361 
362  g.set(0) = cf0[0] ;
363  for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;
364  g.set(nm1s2) = cf0[ n3c*nm1 ] ;
365 
366 // Transformation de Fourier inverse de G
367 //---------------------------------------
368 
369 // FFT inverse
370  fftw_execute(p) ;
371 
372 // Valeurs de f deduites de celles de G
373 //-------------------------------------
374 
375  for ( i = 1; i < nm1s2 ; i++ ) {
376 // ... indice du pt symetrique de psi par rapport a pi/2:
377  int isym = nm1 - i ;
378 
379  double fp = 0.5 * ( g(i) + g(isym) ) ;
380  double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
381  ff0[ n3f*i ] = fp + fm ;
382  ff0[ n3f*isym ] = fp - fm ;
383  }
384 
385 //... cas particuliers:
386  ff0[0] = g(0) + fmoins0 ;
387  ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
388  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
389 
390  } // fin de la boucle sur r
391 
392 //-----------------------------------------------------------------------
393 // partie sin(m phi) avec m impair : transformation cos(l theta) inverse
394 //-----------------------------------------------------------------------
395 
396  j++ ;
397 
398  if ( (j != 1) && (j != n1f-1 ) ) {
399 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
400 // pas nuls
401 
402  for (k=0; k<n3c; k++) {
403 
404  int i0 = n2n3c * j + k ; // indice de depart
405  double* cf0 = cf + i0 ; // tableau des donnees a transformer
406 
407  i0 = n2n3f * j + k ; // indice de depart
408  double* ff0 = ff + i0 ; // tableau resultat
409 
410 // Coefficients impairs de G
411 //--------------------------
412 
413  double c1 = cf0[n3c] ;
414 
415  double som = 0;
416  ff0[n3f] = 0 ;
417  for ( i = 3; i < nt; i += 2 ) {
418  ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
419  som += ff0[ n3f*i ] ;
420  }
421 
422 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
423  double fmoins0 = nm1s2 * c1 + som ;
424 
425 // Coef. impairs de G
426 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
427 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
428  for ( i = 3; i < nt; i += 2 ) {
429  g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
430  }
431 
432 
433 // Coefficients pairs de G
434 //------------------------
435 // Ces coefficients sont egaux aux coefficients pairs du developpement de
436 // f.
437 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
438 // donnait exactement les coef. des cosinus, ce facteur serait 1.
439 
440  g.set(0) = cf0[0] ;
441  for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;
442  g.set(nm1s2) = cf0[ n3c*nm1 ] ;
443 
444 // Transformation de Fourier inverse de G
445 //---------------------------------------
446 
447 // FFT inverse
448  fftw_execute(p) ;
449 
450 // Valeurs de f deduites de celles de G
451 //-------------------------------------
452 
453  for ( i = 1; i < nm1s2 ; i++ ) {
454 // ... indice du pt symetrique de psi par rapport a pi/2:
455  int isym = nm1 - i ;
456 
457  double fp = 0.5 * ( g(i) + g(isym) ) ;
458  double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
459  ff0[ n3f*i ] = fp + fm ;
460  ff0[ n3f*isym ] = fp - fm ;
461  }
462 
463 //... cas particuliers:
464  ff0[0] = g(0) + fmoins0 ;
465  ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
466  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
467  } // fin de la boucle sur r
468 
469  } // fin du cas ou le calcul etait necessaire (i.e. ou les
470  // coef en phi n'etaient pas nuls)
471 
472 // On passe au cas m impair suivant:
473  j+=3 ;
474 
475  } // fin de la boucle sur les cas m impair
476 
477 }
478 }
Lorene prototypes.
Definition: app_hor.h:67