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