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