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