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