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