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