LORENE
FFTW3/citcossinsp.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 citcossinsp_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinsp.C,v 1.4 2014/10/13 08:53:20 j_novak Exp $" ;
24 
25 /*
26  * Transformation inverse sin(2*l*theta) ou cos((2*l+1)*theta) (suivant la
27  * parite de l'indice m en phi) sur le deuxieme indice (theta)
28  * d'un tableau 3-D representant une fonction antisymetrique 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* dimc : tableau du nombre d'elements de cf dans chacune des trois
39  * dimensions.
40  * On doit avoir dimc[1] >= deg[1] = nt.
41  *
42  * double* cf : tableau des coefficients c_l de la fonction definis
43  * comme suit (a r et phi fixes)
44  *
45  * pour m pair:
46  * f(theta) = som_{l=0}^{nt-2} c_l sin( 2 l theta ) .
47  * pour m impair:
48  * f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) .
49  *
50  * L'espace memoire correspondant a ce
51  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
52  * avoir ete alloue avant l'appel a la routine.
53  * Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
54  * le tableau cf comme suit
55  * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
56  * ou j et k sont les indices correspondant a
57  * phi et r respectivement.
58  *
59  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
60  * dimensions.
61  * On doit avoir dimf[1] >= deg[1] = nt.
62  *
63  * Sortie:
64  * -------
65  * double* ff : tableau des valeurs de la fonction aux nt points de
66  * de collocation
67  *
68  * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
69  *
70  * L'espace memoire correspondant a ce
71  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
72  * avoir ete alloue avant l'appel a la routine.
73  * Les valeurs de la fonction sont stokees
74  * dans le tableau ff comme suit
75  * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
76  * ou j et k sont les indices correspondant a
77  * phi et r respectivement.
78  *
79  * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
80  * seul tableau, qui constitue une entree/sortie.
81  *
82  */
83 
84 /*
85  * $Id: citcossinsp.C,v 1.4 2014/10/13 08:53:20 j_novak Exp $
86  * $Log: citcossinsp.C,v $
87  * Revision 1.4 2014/10/13 08:53:20 j_novak
88  * Lorene classes and functions now belong to the namespace Lorene.
89  *
90  * Revision 1.3 2014/10/06 15:18:50 j_novak
91  * Modified #include directives to use c++ syntax.
92  *
93  * Revision 1.2 2013/04/25 15:46:06 j_novak
94  * Added special treatment in the case np = 1, for type_p = NONSYM.
95  *
96  * Revision 1.1 2004/12/21 17:06:03 j_novak
97  * Added all files for using fftw3.
98  *
99  * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
100  * Suppressed the directive #include <malloc.h> for malloc is defined
101  * in <stdlib.h>
102  *
103  * Revision 1.3 2002/10/16 14:36:54 j_novak
104  * Reorganization of #include instructions of standard C++, in order to
105  * use experimental version 3 of gcc.
106  *
107  * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
108  * Modification of declaration of Fortran 77 prototypes for
109  * a better portability (in particular on IBM AIX systems):
110  * All Fortran subroutine names are now written F77_* and are
111  * defined in the new file C++/Include/proto_f77.h.
112  *
113  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
114  * LORENE
115  *
116  * Revision 2.0 1999/02/22 15:42:03 hyc
117  * *** empty log message ***
118  *
119  *
120  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinsp.C,v 1.4 2014/10/13 08:53:20 j_novak Exp $
121  *
122  */
123 
124 
125 // headers du C
126 #include <cstdlib>
127 #include <fftw3.h>
128 
129 //Lorene prototypes
130 #include "tbl.h"
131 
132 // Prototypage des sous-routines utilisees:
133 namespace Lorene {
134 fftw_plan back_fft(int, Tbl*&) ;
135 double* cheb_ini(const int) ;
136 double* chebimp_ini(const int ) ;
137 //*****************************************************************************
138 
139 void citcossinsp(const int* deg, const int* dimc, double* cf, const int* dimf,
140  double* ff)
141 {
142 
143 int i, j, k ;
144 
145 // Dimensions des tableaux ff et cf :
146  int n1f = dimf[0] ;
147  int n2f = dimf[1] ;
148  int n3f = dimf[2] ;
149  int n1c = dimc[0] ;
150  int n2c = dimc[1] ;
151  int n3c = dimc[2] ;
152 
153 // Nombres de degres de liberte en theta :
154  int nt = deg[1] ;
155 
156 // Tests de dimension:
157  if (nt > n2f) {
158  cout << "citcossinsp: nt > n2f : nt = " << nt << " , n2f = "
159  << n2f << endl ;
160  abort () ;
161  exit(-1) ;
162  }
163  if (nt > n2c) {
164  cout << "citcossinsp: nt > n2c : nt = " << nt << " , n2c = "
165  << n2c << endl ;
166  abort () ;
167  exit(-1) ;
168  }
169  if ( (n1f > 1) && (n1c > n1f)) {
170  cout << "citcossinsp: n1c > n1f : n1c = " << n1c << " , n1f = "
171  << n1f << endl ;
172  abort () ;
173  exit(-1) ;
174  }
175  if (n3c > n3f) {
176  cout << "citcossinsp: n3c > n3f : n3c = " << n3c << " , n3f = "
177  << n3f << endl ;
178  abort () ;
179  exit(-1) ;
180  }
181 
182 // Nombre de points pour la FFT:
183  int nm1 = nt - 1;
184  int nm1s2 = nm1 / 2;
185 
186 // Recherche des tables pour la FFT:
187  Tbl* pg = 0x0 ;
188  fftw_plan p = back_fft(nm1, pg) ;
189  Tbl& g = *pg ;
190  double* t1 = new double[nt] ;
191 
192 // Recherche de la table des sin(psi) :
193  double* sinp = cheb_ini(nt);
194 
195 // Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}):
196  double* x = chebimp_ini(nt) ;
197 
198 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
199 // et 0 a dimf[2])
200 
201  int n2n3f = n2f * n3f ;
202  int n2n3c = n2c * n3c ;
203  int borne_phi = n1f-1 ;
204  if (n1f == 1) borne_phi = 1 ;
205 
206 //=======================================================================
207 // Cas m pair
208 //=======================================================================
209 
210  j = 0 ;
211 
212  while (j < borne_phi) { //le dernier coef en phi n'est pas considere
213  // (car nul)
214 
215 //-----------------------------------------------------------------------
216 // partie cos(m phi) avec m pair : transformation sin((2 l) theta) inverse
217 //-----------------------------------------------------------------------
218 
219  for (k=0; k<n3c; k++) {
220 
221  int i0 = n2n3c * j + k ; // indice de depart
222  double* cf0 = cf + i0 ; // tableau des donnees a transformer
223 
224  i0 = n2n3f * j + k ; // indice de depart
225  double* ff0 = ff + i0 ; // tableau resultat
226 
227 
228 /*
229  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
230  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
231  */
232 
233 
234 // Calcul des coefficients de Fourier de la fonction
235 // G(psi) = F+(psi) + F_(psi) sin(psi)
236 // en fonction des coefficients en sin(2l theta) de f:
237 
238 // Coefficients en sinus de G
239 //---------------------------
240 // Ces coefficients sont egaux aux coefficients pairs du developmt. en
241 // sin(2l theta) de f (le facteur -.5 vient de la normalisation
242 // de fftw: si fftw donnait reellement les coefficients en sinus,
243 // il faudrait le remplacer par un +1) :
244 
245  for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;
246 
247 // Coefficients en cosinus de G
248 //-----------------------------
249 // Ces coefficients se deduisent des coefficients impairs du developmt.
250 // en sin(2l theta) de f (le facteur +.25 vient de la normalisation
251 // de fftw: si fftw donnait reellement les coefficients en cosinus,
252 // il faudrait le remplacer par un +.5)
253 
254  g.set(0) = .5 * cf0[n3c] ;
255  for ( i = 1; i < nm1s2; i++ ) {
256  g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
257  }
258  g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
259 
260 
261 // Transformation de Fourier inverse de G
262 //---------------------------------------
263 
264 // FFT inverse
265  fftw_execute(p) ;
266 
267 // Valeurs de f deduites de celles de G
268 //-------------------------------------
269 
270  for ( i = 1; i < nm1s2 ; i++ ) {
271 // ... indice du pt symetrique de psi par rapport a pi/2:
272  int isym = nm1 - i ;
273 
274  double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
275  double fm = 0.5 * ( g(i) - g(isym) ) ;
276  ff0[ n3f*i ] = fp + fm ;
277  ff0[ n3f*isym ] = fp - fm ;
278  }
279 
280 //... cas particuliers:
281  ff0[0] = 0. ;
282  ff0[ n3f*nm1 ] = -2. * g(0) ;
283  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
284 
285  } // fin de la boucle sur r
286 
287 //-----------------------------------------------------------------------
288 // partie sin(m phi) avec m pair : transformation sin( (2 l) theta) inverse
289 //-----------------------------------------------------------------------
290 
291  j++ ;
292 
293  if ( (j != 1) && (j != borne_phi ) ) {
294 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
295 // pas nuls
296 
297  for (k=0; k<n3c; k++) {
298 
299  int i0 = n2n3c * j + k ; // indice de depart
300  double* cf0 = cf + i0 ; // tableau des donnees a transformer
301 
302  i0 = n2n3f * j + k ; // indice de depart
303  double* ff0 = ff + i0 ; // tableau resultat
304 
305 
306 /*
307  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
308  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
309  */
310 
311 
312 // Calcul des coefficients de Fourier de la fonction
313 // G(psi) = F+(psi) + F_(psi) sin(psi)
314 // en fonction des coefficients en sin(2l theta) de f:
315 
316 // Coefficients en sinus de G
317 //---------------------------
318 // Ces coefficients sont egaux aux coefficients pairs du developmt. en
319 // sin(2l theta) de f (le facteur -.5 vient de la normalisation
320 // de fftw: si fftw donnait reellement les coefficients en sinus,
321 // il faudrait le remplacer par un +1) :
322 
323  for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;
324 
325 // Coefficients en cosinus de G
326 //-----------------------------
327 // Ces coefficients se deduisent des coefficients impairs du developmt.
328 // en sin(2l theta) de f (le facteur +.25 vient de la normalisation
329 // de fftw: si fftw donnait reellement les coefficients en cosinus,
330 // il faudrait le remplacer par un +.5)
331 
332  g.set(0) = .5 * cf0[n3c] ;
333  for ( i = 1; i < nm1s2; i++ ) {
334  g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
335  }
336  g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
337 
338 
339 // Transformation de Fourier inverse de G
340 //---------------------------------------
341 
342 // FFT inverse
343  fftw_execute(p) ;
344 
345 // Valeurs de f deduites de celles de G
346 //-------------------------------------
347 
348  for ( i = 1; i < nm1s2 ; i++ ) {
349 // ... indice du pt symetrique de psi par rapport a pi/2:
350  int isym = nm1 - i ;
351 
352  double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
353  double fm = 0.5 * ( g(i) - g(isym) ) ;
354  ff0[ n3f*i ] = fp + fm ;
355  ff0[ n3f*isym ] = fp - fm ;
356  }
357 
358 //... cas particuliers:
359  ff0[0] = 0. ;
360  ff0[ n3f*nm1 ] = -2. * g(0) ;
361  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
362 
363  } // fin de la boucle sur r
364 
365  } // fin du cas ou le calcul etait necessaire (i.e. ou les
366  // coef en phi n'etaient pas nuls)
367 
368 // On passe au cas m pair suivant:
369  j+=3 ;
370 
371  } // fin de la boucle sur les cas m pair
372 
373 //##
374  if (n1f<=3) { // cas m=0 seulement (symetrie axiale)
375  delete [] t1 ;
376  return ;
377  }
378 
379 //=======================================================================
380 // Cas m impair
381 //=======================================================================
382 
383  j = 2 ;
384 
385  while (j < borne_phi) { //le dernier coef en phi n'est pas considere
386  // (car nul)
387 
388 //--------------------------------------------------------------------------
389 // partie cos(m phi) avec m impair : transformation cos((2 l+1) theta) inv.
390 //--------------------------------------------------------------------------
391 
392  for (k=0; k<n3c; k++) {
393 
394  int i0 = n2n3c * j + k ; // indice de depart
395  double* cf0 = cf + i0 ; // tableau des donnees a transformer
396 
397  i0 = n2n3f * j + k ; // indice de depart
398  double* ff0 = ff + i0 ; // tableau resultat
399 
400 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
401 // h(x) := x f(x) (x=cos(theta)) a partir des coefficients de f
402 // (resultat stoke dans le tableau t1 :
403  t1[0] = .5 * cf0[0] ;
404  for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
405  t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
406 
407 /*
408  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
409  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
410  */
411 
412 // Calcul des coefficients de Fourier de la fonction
413 // G(psi) = F+(psi) + F_(psi) sin(psi)
414 // en fonction des coefficients en cos(2l theta) de f:
415 
416 // Coefficients impairs de G
417 //--------------------------
418 
419  double c1 = t1[1] ;
420 
421  double som = 0;
422  ff0[n3f] = 0 ;
423  for ( i = 3; i < nt; i += 2 ) {
424  ff0[ n3f*i ] = t1[i] - c1 ;
425  som += ff0[ n3f*i ] ;
426  }
427 
428 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
429  double fmoins0 = nm1s2 * c1 + som ;
430 
431 // Coef. impairs de G
432 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
433 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
434  for ( i = 3; i < nt; i += 2 ) {
435  g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
436  }
437 
438 
439 // Coefficients pairs de G
440 //------------------------
441 // Ces coefficients sont egaux aux coefficients pairs du developpement de
442 // f.
443 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
444 // donnait exactement les coef. des cosinus, ce facteur serait 1.
445 
446  g.set(0) = t1[0] ;
447  for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
448  g.set(nm1s2) = t1[nm1] ;
449 
450 // Transformation de Fourier inverse de G
451 //---------------------------------------
452 
453 // FFT inverse
454  fftw_execute(p) ;
455 
456 // Valeurs de f deduites de celles de G
457 //-------------------------------------
458 
459  for ( i = 1; i < nm1s2 ; i++ ) {
460 // ... indice du pt symetrique de psi par rapport a pi/2:
461  int isym = nm1 - i ;
462 
463  double fp = 0.5 * ( g(i) + g(isym) ) ;
464  double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
465  ff0[ n3f*i ] = ( fp + fm ) / x[isym] ;
466  ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
467  }
468 
469 //... cas particuliers:
470  ff0[0] = g(0) + fmoins0 ;
471  ff0[ n3f*nm1 ] = 0 ;
472  ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
473 
474  } // fin de la boucle sur r
475 
476 
477 //--------------------------------------------------------------------------
478 // partie sin(m phi) avec m impair : transformation cos((2 l+1) theta) inv.
479 //--------------------------------------------------------------------------
480 
481  j++ ;
482 
483  if ( j != borne_phi ) {
484 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
485 // pas nuls
486 
487  for (k=0; k<n3c; k++) {
488 
489  int i0 = n2n3c * j + k ; // indice de depart
490  double* cf0 = cf + i0 ; // tableau des donnees a transformer
491 
492  i0 = n2n3f * j + k ; // indice de depart
493  double* ff0 = ff + i0 ; // tableau resultat
494 
495 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
496 // h(x) := x f(x) (x=cos(theta)) a partir des coefficients de f
497 // (resultat stoke dans le tableau t1 :
498  t1[0] = .5 * cf0[0] ;
499  for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
500  t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
501 
502 /*
503  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
504  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
505  */
506 
507 // Calcul des coefficients de Fourier de la fonction
508 // G(psi) = F+(psi) + F_(psi) sin(psi)
509 // en fonction des coefficients en cos(2l theta) de f:
510 
511 // Coefficients impairs de G
512 //--------------------------
513 
514  double c1 = t1[1] ;
515 
516  double som = 0;
517  ff0[n3f] = 0 ;
518  for ( i = 3; i < nt; i += 2 ) {
519  ff0[ n3f*i ] = t1[i] - c1 ;
520  som += ff0[ n3f*i ] ;
521  }
522 
523 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
524  double fmoins0 = nm1s2 * c1 + som ;
525 
526 // Coef. impairs de G
527 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
528 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
529  for ( i = 3; i < nt; i += 2 ) {
530  g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
531  }
532 
533 
534 // Coefficients pairs de G
535 //------------------------
536 // Ces coefficients sont egaux aux coefficients pairs du developpement de
537 // f.
538 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
539 // donnait exactement les coef. des cosinus, ce facteur serait 1.
540 
541  g.set(0) = t1[0] ;
542  for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
543  g.set(nm1s2) = t1[nm1] ;
544 
545 // Transformation de Fourier inverse de G
546 //---------------------------------------
547 
548 // FFT inverse
549  fftw_execute(p) ;
550 
551 // Valeurs de f deduites de celles de G
552 //-------------------------------------
553 
554  for ( i = 1; i < nm1s2 ; i++ ) {
555 // ... indice du pt symetrique de psi par rapport a pi/2:
556  int isym = nm1 - i ;
557 
558  double fp = 0.5 * ( g(i) + g(isym) ) ;
559  double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
560  ff0[ n3f*i ] = ( fp + fm ) / x[isym] ;
561  ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
562  }
563 
564 //... cas particuliers:
565  ff0[0] = g(0) + fmoins0 ;
566  ff0[ n3f*nm1 ] = 0 ;
567  ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
568 
569  } // fin de la boucle sur r
570 
571 
572  } // fin du cas ou le calcul etait necessaire (i.e. ou les
573  // coef en phi n'etaient pas nuls)
574 
575 // On passe au cas m impair suivant:
576  j+=3 ;
577 
578  } // fin de la boucle sur les cas m impair
579 
580  delete [] t1 ;
581 }
582 }
Lorene prototypes.
Definition: app_hor.h:64