LORENE
FFT991/cfrchebi.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 /*
27  * Transformation de Tchebyshev (cas rare) sur le troisieme indice (indice
28  * correspondant a r) d'un tableau 3-D decrivant une fonction impaire.
29  * Utilise la routine FFT Fortran FFT991
30  *
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 r est nr = deg[2] et doit etre de la forme
37  * nr = 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[2] >= deg[2] = nr.
41  * NB: pour dimf[0] = 1 (un seul point en phi), la transformation
42  * est bien effectuee.
43  * pour dimf[0] > 1 (plus d'un point en phi), la
44  * transformation n'est effectuee que pour les indices (en phi)
45  * j != 1 et j != dimf[0]-1 (cf. commentaires sur borne_phi).
46  *
47  *
48  * double* ff : tableau des valeurs de la fonction aux nr points de
49  * de collocation
50  *
51  * x_i = sin( pi/2 i/(nr-1) ) 0 <= i <= nr-1
52  *
53  * Les valeurs de la fonction doivent etre stokees dans le
54  * tableau ff comme suit
55  * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ]
56  * ou j et k sont les indices correspondant a phi et theta
57  * respectivement.
58  * L'espace memoire correspondant a ce pointeur doit etre
59  * dimf[0]*dimf[1]*dimf[2] et doit etre alloue avant l'appel a
60  * la routine.
61  *
62  * int* dimc : tableau du nombre d'elements de cc dans chacune des trois
63  * dimensions.
64  * On doit avoir dimc[2] >= deg[2] = nr.
65  *
66  * Sortie:
67  * -------
68  * double* cf : tableau des nr-1 coefficients c_i de la fonction definis
69  * comme suit (a theta et phi fixes)
70  *
71  * f(x) = som_{i=0}^{nr-2} c_i T_{2i+1}(x) ,
72  *
73  * ou T_{2i+1}(x) designe le polynome de Tchebyshev de degre
74  * 2i+1.
75  * Les coefficients c_i (0 <= i <= nr-2) sont stokes dans
76  * le tableau cf comme suit
77  * c_i = cf[ dim[1]*dim[2] * j + dim[2] * k + i ]
78  * ou j et k sont les indices correspondant a phi et theta
79  * respectivement. Par convention, on pose c[nr-1] = 0.
80  * L'espace memoire correspondant a ce pointeur doit etre
81  * dimc[0]*dimc[1]*dimc[2] et doit etre alloue avant l'appel a
82  * la routine.
83  *
84  * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
85  * seul tableau, qui constitue une entree/sortie.
86  *
87  */
88 
89 /*
90  * $Id: cfrchebi.C,v 1.5 2016/12/05 16:18:03 j_novak Exp $
91  * $Log: cfrchebi.C,v $
92  * Revision 1.5 2016/12/05 16:18:03 j_novak
93  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
94  *
95  * Revision 1.4 2014/10/15 12:48:20 j_novak
96  * Corrected namespace declaration.
97  *
98  * Revision 1.3 2014/10/13 08:53:15 j_novak
99  * Lorene classes and functions now belong to the namespace Lorene.
100  *
101  * Revision 1.2 2014/10/06 15:18:44 j_novak
102  * Modified #include directives to use c++ syntax.
103  *
104  * Revision 1.1 2004/12/21 17:06:01 j_novak
105  * Added all files for using fftw3.
106  *
107  * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
108  * Suppressed the directive #include <malloc.h> for malloc is defined
109  * in <stdlib.h>
110  *
111  * Revision 1.3 2002/10/16 14:36:44 j_novak
112  * Reorganization of #include instructions of standard C++, in order to
113  * use experimental version 3 of gcc.
114  *
115  * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
116  * Modification of declaration of Fortran 77 prototypes for
117  * a better portability (in particular on IBM AIX systems):
118  * All Fortran subroutine names are now written F77_* and are
119  * defined in the new file C++/Include/proto_f77.h.
120  *
121  * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
122  * LORENE
123  *
124  * Revision 2.0 1999/02/22 15:48:41 hyc
125  * *** empty log message ***
126  *
127  *
128  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cfrchebi.C,v 1.5 2016/12/05 16:18:03 j_novak Exp $
129  *
130  */
131 
132 
133 // headers du C
134 #include <cassert>
135 #include <cstdlib>
136 
137 // Prototypes of F77 subroutines
138 #include "headcpp.h"
139 #include "proto_f77.h"
140 
141 // Prototypage des sous-routines utilisees:
142 namespace Lorene {
143 int* facto_ini(int ) ;
144 double* trigo_ini(int ) ;
145 double* cheb_ini(const int) ;
146 double* chebimp_ini(const int ) ;
147 
148 //*****************************************************************************
149 
150 void cfrchebi(const int* deg, const int* dimf, double* ff, const int* dimc,
151  double* cf)
152 
153 {
154 
155 int i, j, k ;
156 
157 // Dimensions des tableaux ff et cf :
158  int n1f = dimf[0] ;
159  int n2f = dimf[1] ;
160  int n3f = dimf[2] ;
161  int n1c = dimc[0] ;
162  int n2c = dimc[1] ;
163  int n3c = dimc[2] ;
164 
165 // Nombres de degres de liberte en theta et r :
166  int nr = deg[2] ;
167 
168 // Tests de dimension:
169  if (nr > n3f) {
170  cout << "cfrchebi: nr > n3f : nr = " << nr << " , n3f = "
171  << n3f << endl ;
172  abort () ;
173  exit(-1) ;
174  }
175  if (nr > n3c) {
176  cout << "cfrchebi: nr > n3c : nr = " << nr << " , n3c = "
177  << n3c << endl ;
178  abort () ;
179  exit(-1) ;
180  }
181  if (n1f > n1c) {
182  cout << "cfrchebi: n1f > n1c : n1f = " << n1f << " , n1c = "
183  << n1c << endl ;
184  abort () ;
185  exit(-1) ;
186  }
187  if (n2f > n2c) {
188  cout << "cfrchebi: n2f > n2c : n2f = " << n2f << " , n2c = "
189  << n2c << endl ;
190  abort () ;
191  exit(-1) ;
192  }
193 
194 // Nombre de points pour la FFT:
195  int nm1 = nr - 1;
196  int nm1s2 = nm1 / 2;
197 
198 // Recherche des tables pour la FFT:
199  int* facto = facto_ini(nm1) ;
200  double* trigo = trigo_ini(nm1) ;
201 
202 // Recherche de la table des sin(psi) :
203  double* sinp = cheb_ini(nr);
204 
205 // Recherche de la table des points de collocations x_k :
206  double* x = chebimp_ini(nr);
207 
208 // tableau de travail G et t1
209 // (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
210  double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) );
211  double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
212 
213 // Parametres pour la routine FFT991
214  int jump = 1 ;
215  int inc = 1 ;
216  int lot = 1 ;
217  int isign = - 1 ;
218 
219 // boucle sur phi et theta
220 
221  int n2n3f = n2f * n3f ;
222  int n2n3c = n2c * n3c ;
223 
224 /*
225  * Borne de la boucle sur phi:
226  * si n1f = 1, on effectue la boucle une fois seulement.
227  * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
228  * j=n1f-1 et j=0 ne sont pas consideres car nuls).
229  */
230  int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
231 
232  for (j=0; j< borne_phi; j++) {
233 
234  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
235 
236  for (k=0; k<n2f; k++) {
237 
238  int i0 = n2n3f * j + n3f * k ; // indice de depart
239  double* ff0 = ff + i0 ; // tableau des donnees a transformer
240 
241  i0 = n2n3c * j + n3c * k ; // indice de depart
242  double* cf0 = cf + i0 ; // tableau resultat
243 
244 // Multiplication de la fonction par x (pour la rendre paire)
245 // NB: dans les commentaires qui suivent, on note h(x) = x f(x).
246 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
247 // tableau cf0).
248  cf0[0] = 0 ;
249  for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
250 
251 /*
252  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
253  * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)).
254  */
255 
256 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
257  double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
258 
259 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
260 //---------------------------------------------
261  for ( i = 1; i < nm1s2 ; i++ ) {
262 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
263  int isym = nm1 - i ;
264 // ... indice (dans le tableau cf0) du point x correspondant a psi
265  int ix = nm1 - i ;
266 // ... indice (dans le tableau cf0) du point x correspondant a sym(psi)
267  int ixsym = nm1 - isym ;
268 
269 // ... F+(psi)
270  double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
271 // ... F_(psi) sin(psi)
272  double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
273  g[i] = fp + fms ;
274  g[isym] = fp - fms ;
275  }
276 //... cas particuliers:
277  g[0] = 0.5 * ( cf0[nm1] + cf0[0] );
278  g[nm1s2] = cf0[nm1s2];
279 
280 // Developpement de G en series de Fourier par une FFT
281 //----------------------------------------------------
282 
283  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
284 
285 // Coefficients pairs du developmt. de Tchebyshev de h
286 //----------------------------------------------------
287 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
288 // de G en series de Fourier (le facteur 2. vient de la normalisation
289 // de fft991; si fft991 donnait reellement les coef. en cosinus, il faudrait le
290 // remplacer par un +1.) :
291 
292  cf0[0] = g[0] ;
293  for (i=2; i<nm1; i += 2 ) cf0[i] = 2. * g[i] ;
294  cf0[nm1] = g[nm1] ;
295 
296 // Coefficients impairs du developmt. de Tchebyshev de h
297 //------------------------------------------------------
298 // 1. Coef. c'_k (recurrence amorcee a partir de zero)
299 // Le +4. en facteur de g[i] est du a la normalisation de fft991
300 // (si fft991 donnait reellement les coef. en sinus, il faudrait le
301 // remplacer par un -2.)
302  cf0[1] = 0 ;
303  double som = 0;
304  for ( i = 3; i < nr; i += 2 ) {
305  cf0[i] = cf0[i-2] + 4. * g[i] ;
306  som += cf0[i] ;
307  }
308 
309 // 2. Calcul de c_1 :
310  double c1 = ( fmoins0 - som ) / nm1s2 ;
311 
312 // 3. Coef. c_k avec k impair:
313  cf0[1] = c1 ;
314  for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
315 
316 // Coefficients de f en fonction de ceux de h
317 //-------------------------------------------
318 
319  cf0[0] = 2* cf0[0] ;
320  for (i=1; i<nm1; i++) {
321  cf0[i] = 2 * cf0[i] - cf0[i-1] ;
322  }
323  cf0[nm1] = 0 ;
324 
325 
326  } // fin de la boucle sur theta
327  } // fin de la boucle sur phi
328 
329  // Menage
330  free (t1) ;
331  free (g) ;
332 
333 }
334 }
Lorene prototypes.
Definition: app_hor.h:67