LORENE
FFT991/cftcossinc.C
1 /*
2  * Copyright (c) 1999-2002 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 cftcossinc_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cftcossinc.C,v 1.4 2014/10/15 12:48:20 j_novak Exp $" ;
24 
25 /*
26  * Transformation en sin(l*theta) ou cos(l*theta) (suivant la parite
27  * de l'indice m en phi) sur le deuxieme indice (theta)
28  * d'un tableau 3-D representant une fonction sans symetrie 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* dimf : tableau du nombre d'elements de ff dans chacune des trois
39  * dimensions.
40  * On doit avoir dimf[1] >= deg[1] = nt.
41  *
42  * double* ff : tableau des valeurs de la fonction aux nt points de
43  * de collocation
44  *
45  * theta_l = pi l/(nt-1) 0 <= l <= nt-1
46  *
47  * L'espace memoire correspondant a ce
48  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
49  * etre alloue avant l'appel a la routine.
50  * Les valeurs de la fonction doivent etre stokees
51  * dans le tableau ff comme suit
52  * f( theta_l ) = ff[ dimf[1]*dimf[2] * m + k + dimf[2] * l ]
53  * ou m et k sont les indices correspondant a
54  * phi et r respectivement.
55  * NB: cette routine suppose que la transformation en phi a deja ete
56  * effectuee: ainsi m est un indice de Fourier, non un indice de
57  * point de collocation en phi.
58  *
59  * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
60  * dimensions.
61  * On doit avoir dimc[1] >= deg[1] = nt.
62  * Sortie:
63  * -------
64  * double* cf : tableau des coefficients c_l de la fonction definis
65  * comme suit (a r et phi fixes)
66  *
67  * pour m pair:
68  * f(theta) = som_{l=0}^{nt-1} c_l cos(l theta ) .
69  * pour m impair:
70  * f(theta) = som_{l=0}^{nt-1} c_l sin( l theta ) .
71  *
72  * L'espace memoire correspondant a ce
73  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
74  * etre alloue avant l'appel a la routine.
75  * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
76  * le tableau cf comme suit
77  * c_l = cf[ dimc[1]*dimc[2] * m + k + dimc[2] * l ]
78  * ou m et k sont les indices correspondant a
79  * phi et r respectivement.
80  * Pour m pair, c_0 = c_{nt-1} = 0.
81  * Pour m impair, c_{nt-1} = 0.
82  *
83  * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
84  * seul tableau, qui constitue une entree/sortie.
85  *
86  */
87 
88 /*
89  * $Id: cftcossinc.C,v 1.4 2014/10/15 12:48:20 j_novak Exp $
90  * $Log: cftcossinc.C,v $
91  * Revision 1.4 2014/10/15 12:48:20 j_novak
92  * Corrected namespace declaration.
93  *
94  * Revision 1.3 2014/10/13 08:53:15 j_novak
95  * Lorene classes and functions now belong to the namespace Lorene.
96  *
97  * Revision 1.2 2014/10/06 15:18:45 j_novak
98  * Modified #include directives to use c++ syntax.
99  *
100  * Revision 1.1 2004/12/21 17:06:01 j_novak
101  * Added all files for using fftw3.
102  *
103  * Revision 1.1 2004/11/23 15:13:50 m_forot
104  * Added the bases for the cases without any equatorial symmetry
105  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
106  *
107  *
108  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cftcossinc.C,v 1.4 2014/10/15 12:48:20 j_novak Exp $
109  *
110  */
111 
112 // headers du C
113 #include <cassert>
114 #include <cstdlib>
115 
116 // Prototypes of F77 subroutines
117 #include "headcpp.h"
118 #include "proto_f77.h"
119 
120 // Prototypage des sous-routines utilisees:
121 namespace Lorene {
122 int* facto_ini(int ) ;
123 double* trigo_ini(int ) ;
124 double* cheb_ini(const int) ;
125 double* chebimp_ini(const int ) ;
126 //*****************************************************************************
127 
128 void cftcossinc(const int* deg, const int* dimf, double* ff, const int* dimc,
129  double* cf)
130 {
131 
132 int i, j, k ;
133 
134 // Dimensions des tableaux ff et cf :
135  int n1f = dimf[0] ;
136  int n2f = dimf[1] ;
137  int n3f = dimf[2] ;
138  int n1c = dimc[0] ;
139  int n2c = dimc[1] ;
140  int n3c = dimc[2] ;
141 
142 // Nombre de degres de liberte en theta :
143  int nt = deg[1] ;
144 
145 // Tests de dimension:
146  if (nt > n2f) {
147  cout << "cftcossinc: nt > n2f : nt = " << nt << " , n2f = "
148  << n2f << endl ;
149  abort () ;
150  exit(-1) ;
151  }
152  if (nt > n2c) {
153  cout << "cftcossinc: nt > n2c : nt = " << nt << " , n2c = "
154  << n2c << endl ;
155  abort () ;
156  exit(-1) ;
157  }
158  if (n1f > n1c) {
159  cout << "cftcossinc: n1f > n1c : n1f = " << n1f << " , n1c = "
160  << n1c << endl ;
161  abort () ;
162  exit(-1) ;
163  }
164  if (n3f > n3c) {
165  cout << "cftcossinc: n3f > n3c : n3f = " << n3f << " , n3c = "
166  << n3c << endl ;
167  abort () ;
168  exit(-1) ;
169  }
170 
171 // Nombre de points pour la FFT:
172  int nm1 = nt - 1;
173  int nm1s2 = nm1 / 2;
174 
175 // Recherche des tables pour la FFT:
176  int* facto = facto_ini(nm1) ;
177  double* trigo = trigo_ini(nm1) ;
178 
179 // Recherche de la table des sin(psi) :
180  double* sinp = cheb_ini(nt);
181 
182 // Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}) :
183  double* x = chebimp_ini(nt);
184 
185  // tableau de travail G et t1
186  // (la dimension nm1+2 = nt+1 est exigee par la routine fft991)
187  double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) );
188  double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
189 
190 // Parametres pour la routine FFT991
191  int jump = 1 ;
192  int inc = 1 ;
193  int lot = 1 ;
194  int isign = - 1 ;
195 
196 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
197 // et 0 a dimf[2])
198 
199  int n2n3f = n2f * n3f ;
200  int n2n3c = n2c * n3c ;
201 
202 //=======================================================================
203 // Cas m pair
204 //=======================================================================
205 
206  j = 0 ;
207 
208  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
209  // (car nul)
210 
211 //------------------------------------------------------------------------
212 // partie cos(m phi) avec m pair : transformation en cos(l theta)
213 //------------------------------------------------------------------------
214 
215 
216  for (k=0; k<n3f; k++) {
217 
218  int i0 = n2n3f * j + k ; // indice de depart
219  double* ff0 = ff + i0 ; // tableau des donnees a transformer
220 
221  i0 = n2n3c * j + k ; // indice de depart
222  double* cf0 = cf + i0 ; // tableau resultat
223 
224 
225 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
226  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
227 
228 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
229 //---------------------------------------------
230  for ( i = 1; i < nm1s2 ; i++ ) {
231  int isym = nm1 - i ;
232  int ix = n3f * i ;
233  int ixsym = n3f * isym ;
234 // ... F+(theta)
235  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
236 // ... F_(theta) sin(psi)
237  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
238  g[i] = fp + fms ;
239  g[isym] = fp - fms ;
240  }
241 //... cas particuliers:
242  g[0] = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
243  g[nm1s2] = ff0[ n3f*nm1s2 ];
244 
245 // Developpement de G en series de Fourier par une FFT
246 //----------------------------------------------------
247 
248  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
249 
250 // Coefficients pairs du developmt. cos(l theta) de f
251 //----------------------------------------------------
252 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
253 // de G en series de Fourier (le facteur 2 vient de la normalisation
254 // de fft991) :
255 
256  cf0[0] = g[0] ;
257  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = 2.* g[i] ;
258  cf0[n3c*nm1] = g[nm1] ;
259 
260 // Coefficients impairs du developmt. en cos(l theta) de f
261 //---------------------------------------------------------
262 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
263 // Le +4. en facteur de g[i] est du a la normalisation de fft991
264 // (si fft991 donnait reellement les coef. en sinus, il faudrait le
265 // remplacer par un -2.)
266  cf0[n3c] = 0 ;
267  double som = 0;
268  for ( i = 3; i < nt; i += 2 ) {
269  cf0[n3c*i] = cf0[n3c*(i-2)] + 4. * g[i] ;
270  som += cf0[n3c*i] ;
271  }
272 
273 // 2. Calcul de c_1 :
274  double c1 = ( fmoins0 - som ) / nm1s2 ;
275 
276 // 3. Coef. c_k avec k impair:
277  cf0[n3c] = c1 ;
278  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
279 
280 
281  } // fin de la boucle sur r
282 
283 //--------------------------------------------------------------------
284 // partie sin(m phi) avec m pair : transformation en cos(l theta)
285 //--------------------------------------------------------------------
286 
287  j++ ;
288 
289  if ( (j != 1) && (j != n1f-1 ) ) {
290 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
291 // pas nuls
292 
293  for (k=0; k<n3f; k++) {
294 
295  int i0 = n2n3f * j + k ; // indice de depart
296  double* ff0 = ff + i0 ; // tableau des donnees a transformer
297 
298  i0 = n2n3c * j + k ; // indice de depart
299  double* cf0 = cf + i0 ; // tableau resultat
300 
301 
302 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
303  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
304 
305 // Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
306 //---------------------------------------------
307  for ( i = 1; i < nm1s2 ; i++ ) {
308  int isym = nm1 - i ;
309  int ix = n3f * i ;
310  int ixsym = n3f * isym ;
311 // ... F+(theta)
312  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
313 // ... F_(theta) sin(psi)
314  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
315  g[i] = fp + fms ;
316  g[isym] = fp - fms ;
317  }
318 //... cas particuliers:
319  g[0] = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
320  g[nm1s2] = ff0[ n3f*nm1s2 ];
321 
322 // Developpement de G en series de Fourier par une FFT
323 //----------------------------------------------------
324 
325  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
326 
327 // Coefficients pairs du developmt. cos(l theta) de f
328 //----------------------------------------------------
329 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
330 // de G en series de Fourier (le facteur 2 vient de la normalisation
331 // de fft991) :
332 
333  cf0[0] = g[0] ;
334  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = 2.* g[i] ;
335  cf0[n3c*nm1] = g[nm1] ;
336 
337 // Coefficients impairs du developmt. en cos(l theta) de f
338 //---------------------------------------------------------
339 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
340 // Le +4. en facteur de g[i] est du a la normalisation de fft991
341 // (si fft991 donnait reellement les coef. en sinus, il faudrait le
342 // remplacer par un -2.)
343  cf0[n3c] = 0 ;
344  double som = 0;
345  for ( i = 3; i < nt; i += 2 ) {
346  cf0[n3c*i] = cf0[n3c*(i-2)] + 4. * g[i] ;
347  som += cf0[n3c*i] ;
348  }
349 
350 // 2. Calcul de c_1 :
351  double c1 = ( fmoins0 - som ) / nm1s2 ;
352 
353 // 3. Coef. c_k avec k impair:
354  cf0[n3c] = c1 ;
355  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
356 
357 
358  } // fin de la boucle sur r
359  } // fin du cas ou le calcul etait necessaire (i.e. ou les
360  // coef en phi n'etaient pas nuls)
361 
362 // On passe au cas m pair suivant:
363  j+=3 ;
364 
365  } // fin de la boucle sur les cas m pair
366 
367  if (n1f<=3) { // cas m=0 seulement (symetrie axiale)
368  free (t1) ;
369  free (g) ;
370  return ;
371  }
372 
373 //=======================================================================
374 // Cas m impair
375 //=======================================================================
376 
377  j = 2 ;
378 
379  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
380  // (car nul)
381 
382 //--------------------------------------------------------------------
383 // partie cos(m phi) avec m impair : transformation en sin(l) theta)
384 //--------------------------------------------------------------------
385 
386  for (k=0; k<n3f; k++) {
387 
388  int i0 = n2n3f * j + k ; // indice de depart
389  double* ff0 = ff + i0 ; // tableau des donnees a transformer
390 
391  i0 = n2n3c * j + k ; // indice de depart
392  double* cf0 = cf + i0 ; // tableau resultat
393 
394 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
395  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
396 
397 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
398 //---------------------------------------------
399 
400  for ( i = 1; i < nm1s2 ; i++ ) {
401  int isym = nm1 - i ;
402  int ix = n3f * i ;
403  int ixsym = n3f * isym ;
404 // ... F+(theta)
405  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
406 // ... F_(theta) sin(theta)
407  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
408  g[i] = fp + fms ;
409  g[isym] = fp - fms ;
410  }
411 //... cas particuliers:
412  g[0] = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
413  g[nm1s2] = ff0[ n3f*nm1s2 ];
414 
415 // Developpement de G en series de Fourier par une FFT
416 //----------------------------------------------------
417 
418  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
419 
420 // Coefficients pairs du developmt. sin(l theta) de f
421 //----------------------------------------------------
422 // Ces coefficients sont egaux aux coefficients en sinus du developpement
423 // de G en series de Fourier (le facteur -2 vient de la normalisation
424 // de fft991) :
425 
426  cf0[0] = 0. ;
427  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = -2.* g[i+1] ;
428  cf0[n3c*nm1] = 0. ;
429 
430 // Coefficients impairs du developmt. en sin(l theta) de f
431 //---------------------------------------------------------
432 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
433 // Le +4. en facteur de g[i] est du a la normalisation de fft991
434 // (si fft991 donnait reellement les coef. en sinus, il faudrait le
435 // remplacer par un -2.)
436 
437  cf0[n3c] = 2.* g[0];
438  for ( i = 3; i < nt; i += 2 ) {
439  cf0[n3c*i] = cf0[n3c*(i-2)] + 4. * g[i-1] ;
440  }
441 
442  } // fin de la boucle sur r
443 
444 //------------------------------------------------------------------------
445 // partie sin(m phi) avec m impair : transformation en sin(l theta)
446 //------------------------------------------------------------------------
447 
448  j++ ;
449 
450  if ( j != n1f-1 ) {
451 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
452 // pas nuls
453 
454  for (k=0; k<n3f; k++) {
455 
456  int i0 = n2n3f * j + k ; // indice de depart
457  double* ff0 = ff + i0 ; // tableau des donnees a transformer
458 
459  i0 = n2n3c * j + k ; // indice de depart
460  double* cf0 = cf + i0 ; // tableau resultat
461 
462 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
463  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
464 
465 // Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
466 //---------------------------------------------
467  for ( i = 1; i < nm1s2 ; i++ ) {
468  int isym = nm1 - i ;
469  int ix = n3f * i ;
470  int ixsym = n3f * isym ;
471 // ... F+(theta)
472  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
473 // ... F_(theta) sin(theta)
474  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
475  g[i] = fp + fms ;
476  g[isym] = fp - fms ;
477  }
478 //... cas particuliers:
479  g[0] = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
480  g[nm1s2] = ff0[ n3f*nm1s2 ];
481 
482 // Developpement de G en series de Fourier par une FFT
483 //----------------------------------------------------
484 
485  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
486 
487 // Coefficients pairs du developmt. sin(l theta) de f
488 //----------------------------------------------------
489 // Ces coefficients sont egaux aux coefficients en sinus du developpement
490 // de G en series de Fourier (le facteur -2 vient de la normalisation
491 // de fft991) :
492 
493  cf0[0] = 0. ;
494  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = -2.* g[i+1] ;
495  cf0[n3c*nm1] = 0. ;
496 
497 // Coefficients impairs du developmt. en sin(l theta) de f
498 //---------------------------------------------------------
499 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
500 // Le +4. en facteur de g[i] est du a la normalisation de fft991
501 // (si fft991 donnait reellement les coef. en sinus, il faudrait le
502 // remplacer par un -2.)
503 
504  cf0[n3c] = 2.* g[0];
505  for ( i = 3; i < nt; i += 2 ) {
506  cf0[n3c*i] = cf0[n3c*(i-2)] + 4. * g[i-1] ;
507  }
508 
509  } // fin de la boucle sur r
510 
511  } // fin du cas ou le calcul etait necessaire (i.e. ou les
512  // coef en phi n'etaient pas nuls)
513 
514 
515 // On passe au cas m impair suivant:
516  j+=3 ;
517 
518  } // fin de la boucle sur les cas m impair
519 
520  // Menage
521  free (t1) ;
522  free (g) ;
523 
524 }
525 }
526 
Lorene
Lorene prototypes.
Definition: app_hor.h:64