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