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