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