LORENE
FFTW3/circhebp.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 rare) sur le troisieme indice
28  * (indice correspondant a r) d'un tableau 3-D decrivant une fonction paire.
29  * Utilise 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_{2i}(x) ,
50  *
51  * ou T_{2i}(x) designe le polynome de Tchebyshev de degre 2i.
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 a ce pointeur doit etre
58  * dimc[0]*dimc[1]*dimc[2] et doit etre alloue avant l'appel a
59  * 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 = sin( pi/2 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 avoir ete alloue avant
79  * l'appel a 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  * $Id: circhebp.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
87  * $Log: circhebp.C,v $
88  * Revision 1.4 2016/12/05 16:18:05 j_novak
89  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
90  *
91  * Revision 1.3 2014/10/13 08:53:20 j_novak
92  * Lorene classes and functions now belong to the namespace Lorene.
93  *
94  * Revision 1.2 2014/10/06 15:18:49 j_novak
95  * Modified #include directives to use c++ syntax.
96  *
97  * Revision 1.1 2004/12/21 17:06:02 j_novak
98  * Added all files for using fftw3.
99  *
100  * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
101  * Suppressed the directive #include <malloc.h> for malloc is defined
102  * in <stdlib.h>
103  *
104  * Revision 1.3 2002/10/16 14:36:53 j_novak
105  * Reorganization of #include instructions of standard C++, in order to
106  * use experimental version 3 of gcc.
107  *
108  * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
109  * Modification of declaration of Fortran 77 prototypes for
110  * a better portability (in particular on IBM AIX systems):
111  * All Fortran subroutine names are now written F77_* and are
112  * defined in the new file C++/Include/proto_f77.h.
113  *
114  * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
115  * LORENE
116  *
117  * Revision 2.0 1999/02/22 15:43:29 hyc
118  * *** empty log message ***
119  *
120  *
121  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/circhebp.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
122  *
123  */
124 
125 // headers du C
126 #include <cstdlib>
127 #include <fftw3.h>
128 
129 //Lorene prototypes
130 #include "tbl.h"
131 
132 // Prototypage des sous-routines utilisees:
133 namespace Lorene {
134 fftw_plan back_fft(int, Tbl*&) ;
135 double* cheb_ini(const int) ;
136 //*****************************************************************************
137 
138 void circhebp(const int* deg, const int* dimc, double* cf,
139  const int* dimf, double* ff)
140 
141 {
142 int i, j, k ;
143 
144 // Dimensions des tableaux ff et cf :
145  int n1f = dimf[0] ;
146  int n2f = dimf[1] ;
147  int n3f = dimf[2] ;
148  int n1c = dimc[0] ;
149  int n2c = dimc[1] ;
150  int n3c = dimc[2] ;
151 
152 // Nombres de degres de liberte en r :
153  int nr = deg[2] ;
154 
155 // Tests de dimension:
156  if (nr > n3c) {
157  cout << "circhebp: nr > n3c : nr = " << nr << " , n3c = "
158  << n3c << endl ;
159  abort () ;
160  exit(-1) ;
161  }
162  if (nr > n3f) {
163  cout << "circhebp: nr > n3f : nr = " << nr << " , n3f = "
164  << n3f << endl ;
165  abort () ;
166  exit(-1) ;
167  }
168  if (n1c > n1f) {
169  cout << "circhebp: n1c > n1f : n1c = " << n1c << " , n1f = "
170  << n1f << endl ;
171  abort () ;
172  exit(-1) ;
173  }
174  if (n2c > n2f) {
175  cout << "circhebp: n2c > n2f : n2c = " << n2c << " , n2f = "
176  << n2f << endl ;
177  abort () ;
178  exit(-1) ;
179  }
180 
181 // Nombre de points pour la FFT:
182  int nm1 = nr - 1;
183  int nm1s2 = nm1 / 2;
184 
185 // Recherche des tables pour la FFT:
186  Tbl* pg = 0x0 ;
187  fftw_plan p = back_fft(nm1, pg) ;
188  Tbl& g = *pg ;
189 
190 // Recherche de la table des sin(psi) :
191  double* sinp = cheb_ini(nr);
192 
193 // boucle sur phi et theta
194 
195  int n2n3f = n2f * n3f ;
196  int n2n3c = n2c * n3c ;
197 
198 /*
199  * Borne de la boucle sur phi:
200  * si n1c = 1, on effectue la boucle une fois seulement.
201  * si n1c > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
202  * j=n1c-1 et j=0 ne sont pas consideres car nuls).
203  */
204  int borne_phi = ( n1c > 1 ) ? n1c-1 : 1 ;
205 
206  for (j=0; j< borne_phi; j++) {
207 
208  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
209 
210  for (k=0; k<n2c; k++) {
211 
212  int i0 = n2n3c * j + n3c * k ; // indice de depart
213  double* cf0 = cf + i0 ; // tableau des donnees a transformer
214 
215  i0 = n2n3f * j + n3f * k ; // indice de depart
216  double* ff0 = ff + i0 ; // tableau resultat
217 
218 /*
219  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
220  * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)).
221  */
222 
223 // Calcul des coefficients de Fourier de la fonction
224 // G(psi) = F+(psi) + F_(psi) sin(psi)
225 // en fonction des coefficients de Tchebyshev de f:
226 
227 // Coefficients impairs de G
228 //--------------------------
229 
230  double c1 = cf0[1] ;
231 
232  double som = 0;
233  ff0[1] = 0 ;
234  for ( i = 3; i < nr; i += 2 ) {
235  ff0[i] = cf0[i] - c1 ;
236  som += ff0[i] ;
237  }
238 
239 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
240  double fmoins0 = nm1s2 * c1 + som ;
241 
242 // Coef. impairs de G
243 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
244 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
245  for ( i = 3; i < nr; i += 2 ) {
246  g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ;
247  }
248 
249 
250 // Coefficients pairs de G
251 //------------------------
252 // Ces coefficients sont egaux aux coefficients pairs du developpement de
253 // f.
254 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
255 // donnait exactement les coef. des cosinus, ce facteur serait 1.
256 
257  g.set(0) = cf0[0] ;
258  for (i=1; i<nm1s2; i++) g.set(i) = 0.5 * cf0[2*i] ;
259  g.set(nm1s2) = cf0[nm1] ;
260 
261 // Transformation de Fourier inverse de G
262 //---------------------------------------
263 
264 // FFT inverse
265  fftw_execute(p) ;
266 
267 // Valeurs de f deduites de celles de G
268 //-------------------------------------
269 
270  for ( i = 1; i < nm1s2 ; i++ ) {
271 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
272  int isym = nm1 - i ;
273 // ... indice (dans le tableau ff0) du point x correspondant a psi
274  int ix = nm1 - i ;
275 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
276  int ixsym = nm1 - isym ;
277 
278  double fp = .5 * ( g(i) + g(isym) ) ;
279  double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ;
280 
281  ff0[ix] = fp + fm ;
282  ff0[ixsym] = 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 }
294 }
Lorene prototypes.
Definition: app_hor.h:67