LORENE
FFTW3/circheb.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 inverse (cas fin) sur le troisieme indice
28  * (indice correspondant a r) d'un tableau 3-D
29  * par le biais de la bibliotheque fftw
30  *
31  * Entree:
32  * -------
33  * int* deg : tableau du nombre effectif de degres de liberte dans chacune
34  * des 3 dimensions: le nombre de points de collocation
35  * en r est nr = deg[2] et doit etre de la forme
36  * nr = 2*p + 1
37  * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
38  * dimensions.
39  * On doit avoir dimc[2] >= deg[2] = nr.
40  * NB: pour dimc[0] = 1 (un seul point en phi), la transformation
41  * est bien effectuee.
42  * pour dimc[0] > 1 (plus d'un point en phi), la
43  * transformation n'est effectuee que pour les indices (en phi)
44  * j != 1 et j != dimc[0]-1 (cf. commentaires sur borne_phi).
45  *
46  * double* cf : tableau des coefficients c_i de la fonction definis
47  * comme suit (a theta et phi fixes)
48  *
49  * f(x) = som_{i=0}^{nr-1} c_i T_i(x) ,
50  *
51  * ou T_i(x) designe le polynome de Tchebyshev de degre i.
52  * Les coefficients c_i (0 <= i <= nr-1) doivent etre stokes
53  * dans le tableau cf comme suit
54  * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ]
55  * ou j et k sont les indices correspondant a phi et theta
56  * respectivement.
57  * L'espace memoire correspondant au pointeur cf doit etre
58  * dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant
59  * l'appel a la routine.
60  *
61  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
62  * dimensions.
63  * On doit avoir dimf[2] >= deg[2] = nr.
64  *
65  * Sortie:
66  * -------
67  * double* ff : tableau des valeurs de la fonction aux nr points de
68  * de collocation
69  *
70  * x_i = - cos( pi i/(nr-1) ) 0 <= i <= nr-1
71  *
72  * Les valeurs de la fonction sont stokees dans le
73  * tableau ff comme suit
74  * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ]
75  * ou j et k sont les indices correspondant a phi et theta
76  * respectivement.
77  * L'espace memoire correspondant a ce pointeur doit etre
78  * dimf[0]*dimf[1]*dimf[2] et doit etre alloue avant l'appel a
79  * la routine.
80  *
81  * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
82  * seul tableau, qui constitue une entree/sortie.
83  *
84  */
85 
86 /*
87  * $Id: circheb.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
88  * $Log: circheb.C,v $
89  * Revision 1.4 2016/12/05 16:18:05 j_novak
90  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
91  *
92  * Revision 1.3 2014/10/13 08:53:19 j_novak
93  * Lorene classes and functions now belong to the namespace Lorene.
94  *
95  * Revision 1.2 2014/10/06 15:18:49 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.5 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.4 2002/10/16 14:36:53 j_novak
106  * Reorganization of #include instructions of standard C++, in order to
107  * use experimental version 3 of gcc.
108  *
109  * Revision 1.3 2002/09/09 14:04:22 e_gourgoulhon
110  *
111  * Correction of an error : fft991_ -> F77_fft991
112  *
113  * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
114  * Modification of declaration of Fortran 77 prototypes for
115  * a better portability (in particular on IBM AIX systems):
116  * All Fortran subroutine names are now written F77_* and are
117  * defined in the new file C++/Include/proto_f77.h.
118  *
119  * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
120  * LORENE
121  *
122  * Revision 2.0 1999/02/22 15:43:47 hyc
123  * *** empty log message ***
124  *
125  *
126  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/circheb.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
127  *
128  */
129 
130 // headers du C
131 #include <cassert>
132 #include <cstdlib>
133 #include <fftw3.h>
134 
135 //Lorene prototypes
136 #include "tbl.h"
137 
138 // Prototypage des sous-routines utilisees:
139 namespace Lorene {
140 fftw_plan back_fft(int, Tbl*&) ;
141 double* cheb_ini(const int) ;
142 //*****************************************************************************
143 
144 void circheb(const int* deg, const int* dimc, double* cf, const int* dimf,
145  double* ff)
146 
147 {
148 int i, j, k ;
149 
150 // Dimensions des tableaux ff et cf :
151  int n1f = dimf[0] ;
152  int n2f = dimf[1] ;
153  int n3f = dimf[2] ;
154  int n1c = dimc[0] ;
155  int n2c = dimc[1] ;
156  int n3c = dimc[2] ;
157 
158 // Nombres de degres de liberte en r :
159  int nr = deg[2] ;
160 
161 // Tests de dimension:
162  if (nr > n3c) {
163  cout << "circheb: nr > n3c : nr = " << nr << " , n3c = "
164  << n3c << endl ;
165  abort () ;
166  exit(-1) ;
167  }
168  if (nr > n3f) {
169  cout << "circheb: nr > n3f : nr = " << nr << " , n3f = "
170  << n3f << endl ;
171  abort () ;
172  exit(-1) ;
173  }
174  if (n1c > n1f) {
175  cout << "circheb: n1c > n1f : n1c = " << n1c << " , n1f = "
176  << n1f << endl ;
177  abort () ;
178  exit(-1) ;
179  }
180  if (n2c > n2f) {
181  cout << "circheb: n2c > n2f : n2c = " << n2c << " , n2f = "
182  << n2f << endl ;
183  abort () ;
184  exit(-1) ;
185  }
186 
187 // Nombre de points pour la FFT inverse:
188  int nm1 = nr - 1;
189  int nm1s2 = nm1 / 2;
190 
191 // Recherche des tables pour la FFT inverse:
192  Tbl* pg = 0x0 ;
193  fftw_plan p = back_fft(nm1, pg) ;
194  Tbl& g = *pg ;
195 
196 // Recherche de la table des sin(psi) :
197  double* sinp = cheb_ini(nr);
198 
199 // boucle sur phi et theta
200 
201  int n2n3f = n2f * n3f ;
202  int n2n3c = n2c * n3c ;
203 
204 /*
205  * Borne de la boucle sur phi:
206  * si n1c = 1, on effectue la boucle une fois seulement.
207  * si n1c > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
208  * j=n1c-1 et j=0 ne sont pas consideres car nuls).
209  */
210  int borne_phi = ( n1c > 1 ) ? n1c-1 : 1 ;
211 
212  for (j=0; j< borne_phi; j++) {
213 
214  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
215 
216  for (k=0; k<n2c; k++) {
217 
218  int i0 = n2n3c * j + n3c * k ; // indice de depart
219  double* cf0 = cf + i0 ; // tableau des donnees a transformer
220 
221  i0 = n2n3f * j + n3f * k ; // indice de depart
222  double* ff0 = ff + i0 ; // tableau resultat
223 
224 /*
225  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
226  * reliee a x par x = - cos(psi) et F(psi) = f(x(psi)).
227  */
228 
229 // Calcul des coefficients de Fourier de la fonction
230 // G(psi) = F+(psi) + F_(psi) sin(psi)
231 // en fonction des coefficients de Tchebyshev de f:
232 
233 // Coefficients impairs de G
234 //--------------------------
235 
236  double c1 = cf0[1] ;
237 
238  double som = 0;
239  ff0[1] = 0 ;
240  for ( i = 3; i < nr; i += 2 ) {
241  ff0[i] = cf0[i] - c1 ;
242  som += ff0[i] ;
243  }
244 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
245  double fmoins0 = - nm1s2 * c1 - som ;
246 
247 // Coef. impairs de G
248 // NB: le facteur -0.25 est du a la normalisation de fftw; si fftw
249 // donnait exactement les coef. des sinus, ce facteur serait +0.5.
250  for ( i = 3; i < nr; i += 2 ) {
251  g.set(nm1-i/2) = -0.25 * ( ff0[i] - ff0[i-2] ) ;
252  }
253 
254 // Coefficients pairs de G
255 //------------------------
256 // Ces coefficients sont egaux aux coefficients pairs du developpement de
257 // f en polynomes de Tchebyshev.
258 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
259 // donnait exactement les coef. des cosinus, ce facteur serait 1.
260 
261  g.set(0) = cf0[0] ;
262  for (i=1; i<nm1s2; i ++ ) g.set(i) = 0.5 * cf0[2*i] ;
263  g.set(nm1s2) = cf0[nm1] ;
264 
265 // Transformation de Fourier inverse de G
266 //---------------------------------------
267 
268 // FFT inverse
269  fftw_execute(p) ;
270 
271 // Valeurs de f deduites de celles de G
272 //-------------------------------------
273 
274  for ( i = 1; i < nm1s2 ; i++ ) {
275 // ... indice du pt symetrique de psi par rapport a pi/2:
276  int isym = nm1 - i ;
277 
278  double fp = .5 * ( g(i) + g(isym) ) ;
279  double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ;
280 
281  ff0[i] = fp + fm ;
282  ff0[isym] = fp - fm ;
283  }
284 
285 //... cas particuliers:
286  ff0[0] = g(0) + fmoins0 ;
287  ff0[nm1] = g(0) - fmoins0 ;
288  ff0[nm1s2] = g(nm1s2) ;
289 
290  } // fin de la boucle sur theta
291  } // fin de la boucle sur phi
292 }
293 }
Lorene prototypes.
Definition: app_hor.h:67