LORENE
FFT991/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
23char circheb_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/circheb.C,v 1.4 2014/10/15 12:48:21 j_novak Exp $" ;
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 routine FFT Fortran FFT991
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 3^q 5^r + 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 2014/10/15 12:48:21 j_novak Exp $
88 * $Log: circheb.C,v $
89 * Revision 1.4 2014/10/15 12:48:21 j_novak
90 * Corrected namespace declaration.
91 *
92 * Revision 1.3 2014/10/13 08:53:16 j_novak
93 * Lorene classes and functions now belong to the namespace Lorene.
94 *
95 * Revision 1.2 2014/10/06 15:18:45 j_novak
96 * Modified #include directives to use c++ syntax.
97 *
98 * Revision 1.1 2004/12/21 17:06:01 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/FFT991/circheb.C,v 1.4 2014/10/15 12:48:21 j_novak Exp $
127 *
128 */
129
130// headers du C
131#include <cassert>
132#include <cstdlib>
133
134// Prototypes of F77 subroutines
135#include "headcpp.h"
136#include "proto_f77.h"
137
138// Prototypage des sous-routines utilisees:
139namespace Lorene {
140int* facto_ini(int ) ;
141double* trigo_ini(int ) ;
142double* cheb_ini(const int) ;
143//*****************************************************************************
144
145void circheb(const int* deg, const int* dimc, double* cf, const int* dimf,
146 double* ff)
147
148{
149
150int 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 r :
161 int nr = deg[2] ;
162
163// Tests de dimension:
164 if (nr > n3c) {
165 cout << "circheb: nr > n3c : nr = " << nr << " , n3c = "
166 << n3c << endl ;
167 abort () ;
168 exit(-1) ;
169 }
170 if (nr > n3f) {
171 cout << "circheb: nr > n3f : nr = " << nr << " , n3f = "
172 << n3f << endl ;
173 abort () ;
174 exit(-1) ;
175 }
176 if (n1c > n1f) {
177 cout << "circheb: n1c > n1f : n1c = " << n1c << " , n1f = "
178 << n1f << endl ;
179 abort () ;
180 exit(-1) ;
181 }
182 if (n2c > n2f) {
183 cout << "circheb: n2c > n2f : n2c = " << n2c << " , n2f = "
184 << n2f << endl ;
185 abort () ;
186 exit(-1) ;
187 }
188
189// Nombre de points pour la FFT inverse:
190 int nm1 = nr - 1;
191 int nm1s2 = nm1 / 2;
192
193// Recherche des tables pour la FFT inverse:
194 int* facto = facto_ini(nm1) ;
195 double* trigo = trigo_ini(nm1) ;
196
197// Recherche de la table des sin(psi) :
198 double* sinp = cheb_ini(nr);
199
200 // tableau de travail t1 et g
201 // (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
202 double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
203 double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
204
205// Parametres pour la routine FFT991
206 int jump = 1 ;
207 int inc = 1 ;
208 int lot = 1 ;
209 int isign = 1 ;
210
211// boucle sur phi et theta
212
213 int n2n3f = n2f * n3f ;
214 int n2n3c = n2c * n3c ;
215
216/*
217 * Borne de la boucle sur phi:
218 * si n1c = 1, on effectue la boucle une fois seulement.
219 * si n1c > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
220 * j=n1c-1 et j=0 ne sont pas consideres car nuls).
221 */
222 int borne_phi = ( n1c > 1 ) ? n1c-1 : 1 ;
223
224 for (j=0; j< borne_phi; j++) {
225
226 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
227
228 for (k=0; k<n2c; k++) {
229
230 int i0 = n2n3c * j + n3c * k ; // indice de depart
231 double* cf0 = cf + i0 ; // tableau des donnees a transformer
232
233 i0 = n2n3f * j + n3f * k ; // indice de depart
234 double* ff0 = ff + i0 ; // tableau resultat
235
236/*
237 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
238 * reliee a x par x = - cos(psi) et F(psi) = f(x(psi)).
239 */
240
241// Calcul des coefficients de Fourier de la fonction
242// G(psi) = F+(psi) + F_(psi) sin(psi)
243// en fonction des coefficients de Tchebyshev de f:
244
245// Coefficients impairs de G
246//--------------------------
247
248 double c1 = cf0[1] ;
249
250 double som = 0;
251 ff0[1] = 0 ;
252 for ( i = 3; i < nr; i += 2 ) {
253 ff0[i] = cf0[i] - c1 ;
254 som += ff0[i] ;
255 }
256// Valeur en psi=0 de la partie antisymetrique de F, F_ :
257 double fmoins0 = - nm1s2 * c1 - som ;
258
259// Coef. impairs de G
260// NB: le facteur -0.25 est du a la normalisation de fft991; si fft991
261// donnait exactement les coef. des sinus, ce facteur serait +0.5.
262 g[1] = 0 ;
263 for ( i = 3; i < nr; i += 2 ) {
264 g[i] = -0.25 * ( ff0[i] - ff0[i-2] ) ;
265 }
266 g[nr] = 0 ;
267
268// Coefficients pairs de G
269//------------------------
270// Ces coefficients sont egaux aux coefficients pairs du developpement de
271// f en polynomes de Tchebyshev.
272// NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
273// donnait exactement les coef. des cosinus, ce facteur serait 1.
274
275 g[0] = cf0[0] ;
276 for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * cf0[i] ;
277 g[nm1] = cf0[nm1] ;
278
279// Transformation de Fourier inverse de G
280//---------------------------------------
281
282// FFT inverse
283 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
284
285// Valeurs de f deduites de celles de G
286//-------------------------------------
287
288 for ( i = 1; i < nm1s2 ; i++ ) {
289// ... indice du pt symetrique de psi par rapport a pi/2:
290 int isym = nm1 - i ;
291
292 double fp = .5 * ( g[i] + g[isym] ) ;
293 double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
294 ff0[i] = fp + fm ;
295 ff0[isym] = fp - fm ;
296 }
297
298//... cas particuliers:
299 ff0[0] = g[0] + fmoins0 ;
300 ff0[nm1] = g[0] - fmoins0 ;
301 ff0[nm1s2] = g[nm1s2] ;
302
303 } // fin de la boucle sur theta
304 } // fin de la boucle sur phi
305
306 // Menage
307 free (t1) ;
308 free (g) ;
309}
310}
Lorene prototypes.
Definition app_hor.h:64