LORENE
FFT991/circhebpip.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 quelconque.
29 * Utilise 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 * Si l pair f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) ,
50 * Si l impair f(x) = som_{i=0}^{nr-1} c_i T_{2i+1}(x) ,
51 *
52 * ou T_{i}(x) designe le polynome de Tchebyshev de degre i.
53 * Les coefficients c_i (0 <= i <= nr-1) 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: circhebpip.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
88 * $Log: circhebpip.C,v $
89 * Revision 1.4 2016/12/05 16:18:04 j_novak
90 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
91 *
92 * Revision 1.3 2014/10/13 08:53:17 j_novak
93 * Lorene classes and functions now belong to the namespace Lorene.
94 *
95 * Revision 1.2 2014/10/06 15:18:46 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.1 2004/11/23 15:13:50 m_forot
102 * Added the bases for the cases without any equatorial symmetry
103 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
104 *
105 *
106 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/circhebpip.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
107 *
108 */
109
110// headers du C
111#include <cassert>
112#include <cstdlib>
113
114#include "headcpp.h"
115
116// Prototypes of F77 subroutines
117#include "proto_f77.h"
118
119// Prototypage des sous-routines utilisees:
120namespace Lorene {
121int* facto_ini(int ) ;
122double* trigo_ini(int ) ;
123double* cheb_ini(const int) ;
124double* chebimp_ini(const int ) ;
125//*****************************************************************************
126
127void circhebpip(const int* deg, const int* dimc, double* cf,
128 const int* dimf, double* ff)
129
130{
131
132int i, j, k ;
133
134// Dimensions des tableaux ff et cf :
135 int n1f = dimf[0] ;
136 int n2f = dimf[1] ;
137 int n3f = dimf[2] ;
138 int n1c = dimc[0] ;
139 int n2c = dimc[1] ;
140 int n3c = dimc[2] ;
141
142// Nombres de degres de liberte en r :
143 int nr = deg[2] ;
144
145// Tests de dimension:
146 if (nr > n3c) {
147 cout << "circhebpip: nr > n3c : nr = " << nr << " , n3c = "
148 << n3c << endl ;
149 abort () ;
150 exit(-1) ;
151 }
152 if (nr > n3f) {
153 cout << "circhebpip: nr > n3f : nr = " << nr << " , n3f = "
154 << n3f << endl ;
155 abort () ;
156 exit(-1) ;
157 }
158 if (n1c > n1f) {
159 cout << "circhebpip: n1c > n1f : n1c = " << n1c << " , n1f = "
160 << n1f << endl ;
161 abort () ;
162 exit(-1) ;
163 }
164 if (n2c > n2f) {
165 cout << "circhebpip: n2c > n2f : n2c = " << n2c << " , n2f = "
166 << n2f << endl ;
167 abort () ;
168 exit(-1) ;
169 }
170
171// Nombre de points pour la FFT:
172 int nm1 = nr - 1;
173 int nm1s2 = nm1 / 2;
174
175// Recherche des tables pour la FFT:
176 int* facto = facto_ini(nm1) ;
177 double* trigo = trigo_ini(nm1) ;
178
179// Recherche de la table des sin(psi) :
180 double* sinp = cheb_ini(nr);
181
182// Recherche de la table des points de collocations x_k :
183 double* x = chebimp_ini(nr);
184
185 // tableau de travail t1 et g
186 // (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
187 double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
188 double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
189
190// Parametres pour la routine FFT991
191 int jump = 1 ;
192 int inc = 1 ;
193 int lot = 1 ;
194 int isign = 1 ;
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
214 /************ Cas l pair **********/
215
216 for (k=0; k<n2c; k+=2) {
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// Calcul des coefficients de Fourier de la fonction
226// G(psi) = F+(theta) + F_(theta) sin(theta)
227// en fonction des coefficients de Tchebyshev de f:
228
229// Coefficients impairs de G
230//--------------------------
231
232 double c1 = cf0[1] ;
233
234 double som = 0;
235 ff0[1] = 0 ;
236 for ( i = 3; i < nr; i += 2 ) {
237 ff0[i] = cf0[i] - c1 ;
238 som += ff0[i] ;
239 }
240
241// Valeur en theta=0 de la partie antisymetrique de F, F_ :
242 double fmoins0 = nm1s2 * c1 + som ;
243
244// Coef. impairs de G
245// NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
246// donnait exactement les coef. des sinus, ce facteur serait -0.5.
247 g[1] = 0 ;
248 for ( i = 3; i < nr; i += 2 ) {
249 g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
250 }
251 g[nr] = 0 ;
252
253
254// Coefficients pairs de G
255//------------------------
256// Ces coefficients sont egaux aux coefficients pairs du developpement de
257// f.
258// NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
259// donnait exactement les coef. des cosinus, ce facteur serait 1.
260
261 g[0] = cf0[0] ;
262 for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * cf0[i] ;
263 g[nm1] = cf0[nm1] ;
264
265// Transformation de Fourier inverse de G
266//---------------------------------------
267
268// FFT inverse
269 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
270
271// Valeurs de f deduites de celles de G
272//-------------------------------------
273
274 for ( i = 1; i < nm1s2 ; i++ ) {
275// ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
276 int isym = nm1 - i ;
277// ... indice (dans le tableau ff0) du point x correspondant a theta
278 int ix = nm1 - i ;
279// ... indice (dans le tableau ff0) du point x correspondant a sym(theta)
280 int ixsym = nm1 - isym ;
281
282 double fp = .5 * ( g[i] + g[isym] ) ;
283 double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
284
285 ff0[ix] = fp + fm ;
286 ff0[ixsym] = fp - fm ;
287 }
288
289//... cas particuliers:
290 ff0[0] = g[0] - fmoins0 ;
291 ff0[nm1] = g[0] + fmoins0 ;
292 ff0[nm1s2] = g[nm1s2] ;
293
294 } // fin de la boucle sur theta
295
296 /*********** Cas l impair **********/
297
298 for (k=1; k<n2c; k+=2) {
299
300 int i0 = n2n3c * j + n3c * k ; // indice de depart
301 double* cf0 = cf + i0 ; // tableau des donnees a transformer
302
303 i0 = n2n3f * j + n3f * k ; // indice de depart
304 double* ff0 = ff + i0 ; // tableau resultat
305
306// Calcul des coefficients du developpement en T_{2i}(x) de la fonction
307// h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le
308// tableau t1 :
309 t1[0] = .5 * cf0[0] ;
310 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ;
311 t1[nm1] = .5 * cf0[nr-2] ;
312
313
314// Calcul des coefficients de Fourier de la fonction
315// G(psi) = F+(theta) + F_(theta) sin(theta)
316// en fonction des coefficients de Tchebyshev de f:
317
318// Coefficients impairs de G
319//--------------------------
320
321 double c1 = t1[1] ;
322
323 double som = 0;
324 ff0[1] = 0 ;
325 for ( i = 3; i < nr; i += 2 ) {
326 ff0[i] = t1[i] - c1 ;
327 som += ff0[i] ;
328 }
329
330// Valeur en theta=0 de la partie antisymetrique de F, F_ :
331 double fmoins0 = nm1s2 * c1 + som ;
332
333// Coef. impairs de G
334// NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
335// donnait exactement les coef. des sinus, ce facteur serait -0.5.
336 g[1] = 0 ;
337 for ( i = 3; i < nr; i += 2 ) {
338 g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
339 }
340 g[nr] = 0 ;
341
342
343// Coefficients pairs de G
344//------------------------
345// Ces coefficients sont egaux aux coefficients pairs du developpement de
346// f.
347// NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
348// donnait exactement les coef. des cosinus, ce facteur serait 1.
349
350 g[0] = t1[0] ;
351 for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * t1[i] ;
352 g[nm1] = t1[nm1] ;
353
354// Transformation de Fourier inverse de G
355//---------------------------------------
356
357// FFT inverse
358 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
359
360// Valeurs de f deduites de celles de G
361//-------------------------------------
362
363 for ( i = 1; i < nm1s2 ; i++ ) {
364// ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
365 int isym = nm1 - i ;
366// ... indice (dans le tableau ff0) du point x correspondant a theta
367 int ix = nm1 - i ;
368// ... indice (dans le tableau ff0) du point x correspondant a sym(theta)
369 int ixsym = nm1 - isym ;
370
371 double fp = .5 * ( g[i] + g[isym] ) ;
372 double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
373
374 ff0[ix] = ( fp + fm ) / x[ix];
375 ff0[ixsym] = ( fp - fm ) / x[ixsym] ;
376 }
377
378//... cas particuliers:
379 ff0[0] = 0 ;
380 ff0[nm1] = g[0] + fmoins0 ;
381 ff0[nm1s2] = g[nm1s2] / x[nm1s2] ;
382
383 } // fin de la boucle sur theta
384
385
386 } // fin de la boucle sur phi
387
388 // Menage
389 free (t1) ;
390 free (g) ;
391
392}
393}
394
Lorene prototypes.
Definition app_hor.h:67
Coord x
x coordinate centered on the grid
Definition map.h:738
Coord sinp
Definition map.h:735