LORENE
FFTW3/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 bibliotheque fftw.
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 + 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 2016/12/05 16:18:04 j_novak Exp $
91  * $Log: cfrchebi.C,v $
92  * Revision 1.4 2016/12/05 16:18:04 j_novak
93  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
94  *
95  * Revision 1.3 2014/10/13 08:53:18 j_novak
96  * Lorene classes and functions now belong to the namespace Lorene.
97  *
98  * Revision 1.2 2014/10/06 15:18:47 j_novak
99  * Modified #include directives to use c++ syntax.
100  *
101  * Revision 1.1 2004/12/21 17:06:02 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/FFTW3/cfrchebi.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
126  *
127  */
128 
129 
130 // headers du C
131 #include <cstdlib>
132 #include <fftw3.h>
133 
134 //Lorene prototypes
135 #include "tbl.h"
136 
137 // Prototypage des sous-routines utilisees:
138 namespace Lorene {
139 fftw_plan prepare_fft(int, Tbl*&) ;
140 double* cheb_ini(const int) ;
141 double* chebimp_ini(const int ) ;
142 
143 //*****************************************************************************
144 
145 void cfrchebi(const int* deg, const int* dimf, double* ff, const int* dimc,
146  double* cf)
147 
148 {
149 
150 int i, j, k ;
151 
152 // Dimensions des tableaux ff et cf :
153  int n1f = dimf[0] ;
154  int n2f = dimf[1] ;
155  int n3f = dimf[2] ;
156  int n1c = dimc[0] ;
157  int n2c = dimc[1] ;
158  int n3c = dimc[2] ;
159 
160 // Nombres de degres de liberte en theta et r :
161  int nr = deg[2] ;
162 
163 // Tests de dimension:
164  if (nr > n3f) {
165  cout << "cfrchebi: nr > n3f : nr = " << nr << " , n3f = "
166  << n3f << endl ;
167  abort () ;
168  exit(-1) ;
169  }
170  if (nr > n3c) {
171  cout << "cfrchebi: nr > n3c : nr = " << nr << " , n3c = "
172  << n3c << endl ;
173  abort () ;
174  exit(-1) ;
175  }
176  if (n1f > n1c) {
177  cout << "cfrchebi: n1f > n1c : n1f = " << n1f << " , n1c = "
178  << n1c << endl ;
179  abort () ;
180  exit(-1) ;
181  }
182  if (n2f > n2c) {
183  cout << "cfrchebi: n2f > n2c : n2f = " << n2f << " , n2c = "
184  << n2c << endl ;
185  abort () ;
186  exit(-1) ;
187  }
188 
189 // Nombre de points pour la FFT:
190  int nm1 = nr - 1;
191  int nm1s2 = nm1 / 2;
192 
193 // Recherche des tables pour la FFT:
194  Tbl* pg = 0x0 ;
195  fftw_plan p = prepare_fft(nm1, pg) ;
196  Tbl& g = *pg ;
197 
198 // Recherche de la table des sin(psi) :
199  double* sinp = cheb_ini(nr);
200 
201 // Recherche de la table des points de collocations x_k :
202  double* x = chebimp_ini(nr);
203 
204 // boucle sur phi et theta
205 
206  int n2n3f = n2f * n3f ;
207  int n2n3c = n2c * n3c ;
208 
209 /*
210  * Borne de la boucle sur phi:
211  * si n1f = 1, on effectue la boucle une fois seulement.
212  * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
213  * j=n1f-1 et j=0 ne sont pas consideres car nuls).
214  */
215  int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
216 
217  for (j=0; j< borne_phi; j++) {
218 
219  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
220 
221  for (k=0; k<n2f; k++) {
222 
223  int i0 = n2n3f * j + n3f * k ; // indice de depart
224  double* ff0 = ff + i0 ; // tableau des donnees a transformer
225 
226  i0 = n2n3c * j + n3c * k ; // indice de depart
227  double* cf0 = cf + i0 ; // tableau resultat
228 
229 // Multiplication de la fonction par x (pour la rendre paire)
230 // NB: dans les commentaires qui suivent, on note h(x) = x f(x).
231 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
232 // tableau cf0).
233  cf0[0] = 0 ;
234  for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
235 
236 /*
237  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
238  * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)).
239  */
240 
241 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
242  double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
243 
244 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
245 //---------------------------------------------
246  for ( i = 1; i < nm1s2 ; i++ ) {
247 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
248  int isym = nm1 - i ;
249 // ... indice (dans le tableau cf0) du point x correspondant a psi
250  int ix = nm1 - i ;
251 // ... indice (dans le tableau cf0) du point x correspondant a sym(psi)
252  int ixsym = nm1 - isym ;
253 
254 // ... F+(psi)
255  double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
256 // ... F_(psi) sin(psi)
257  double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
258  g.set(i) = fp + fms ;
259  g.set(isym) = fp - fms ;
260  }
261 //... cas particuliers:
262  g.set(0) = 0.5 * ( cf0[nm1] + cf0[0] );
263  g.set(nm1s2) = cf0[nm1s2];
264 
265 // Developpement de G en series de Fourier par une FFT
266 //----------------------------------------------------
267 
268  fftw_execute(p) ;
269 
270 // Coefficients pairs du developmt. de Tchebyshev de h
271 //----------------------------------------------------
272 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
273 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
274 // de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
275 // remplacer par un +1.) :
276 
277  double fac = 2./double(nm1) ;
278  cf0[0] = g(0) / double(nm1) ;
279  for (i=2; i<nm1; i += 2 ) cf0[i] = fac * g(i/2) ;
280  cf0[nm1] = g(nm1s2) / double(nm1) ;
281 
282 // Coefficients impairs du developmt. de Tchebyshev de h
283 //------------------------------------------------------
284 // 1. Coef. c'_k (recurrence amorcee a partir de zero)
285 // Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
286 // (si fftw donnait reellement les coef. en sinus, il faudrait le
287 // remplacer par un -2.)
288  fac *= 2 ;
289  cf0[1] = 0 ;
290  double som = 0;
291  for ( i = 3; i < nr; i += 2 ) {
292  cf0[i] = cf0[i-2] + fac * g(nm1 - i/2) ;
293  som += cf0[i] ;
294  }
295 
296 // 2. Calcul de c_1 :
297  double c1 = ( fmoins0 - som ) / nm1s2 ;
298 
299 // 3. Coef. c_k avec k impair:
300  cf0[1] = c1 ;
301  for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
302 
303 // Coefficients de f en fonction de ceux de h
304 //-------------------------------------------
305 
306  cf0[0] = 2* cf0[0] ;
307  for (i=1; i<nm1; i++) {
308  cf0[i] = 2 * cf0[i] - cf0[i-1] ;
309  }
310  cf0[nm1] = 0 ;
311 
312 
313  } // fin de la boucle sur theta
314  } // fin de la boucle sur phi
315 
316 
317 }
318 }
Lorene prototypes.
Definition: app_hor.h:67