LORENE
FFT991/cfrchebpii.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 (cas rare) sur le troisieme indice (indice
28 * correspondant a r) d'un tableau 3-D decrivant une fonction quelconque.
29 * Utilise la routine FFT Fortran FFT991
30 *
31 *
32 * Entree:
33 * -------
34 * int* deg : tableau du nombre effectif de degres de liberte dans chacune
35 * des 3 dimensions: le nombre de points de collocation
36 * en r est nr = deg[2] et doit etre de la forme
37 * nr = 2^p 3^q 5^r + 1
38 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
39 * dimensions.
40 * On doit avoir dimf[2] >= deg[2] = nr.
41 * NB: pour dimf[0] = 1 (un seul point en phi), la transformation
42 * est bien effectuee.
43 * pour dimf[0] > 1 (plus d'un point en phi), la
44 * transformation n'est effectuee que pour les indices (en phi)
45 * j != 1 et j != dimf[0]-1 (cf. commentaires sur borne_phi).
46 *
47 * double* ff : tableau des valeurs de la fonction aux nr points de
48 * de collocation
49 *
50 * x_i = sin( pi/2 i/(nr-1) ) 0 <= i <= nr-1
51 *
52 * Les valeurs de la fonction doivent etre stokees dans le
53 * tableau ff comme suit
54 * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[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 * dimf[0]*dimf[1]*dimf[2] et doit etre alloue avant l'appel a
59 * la routine.
60 *
61 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
62 * dimensions.
63 * On doit avoir dimc[2] >= deg[2] = nr.
64 *
65 * Sortie:
66 * -------
67 * double* cf : tableau des nr coefficients c_i de la fonction definis
68 * comme suit (a theta et phi fixes)
69 *
70 * Si l impair f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) ,
71 * Si l pair f(x) = som_{i=0}^{nr-1} c_i T_{2i+1}(x)
72 *
73 * ou T_{i}(x) designe le polynome de Tchebyshev de degre i.
74 * Les coefficients c_i (0 <= i <= nr-1) sont stokes dans
75 * le tableau cf comme suit
76 * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ]
77 * ou j et k sont les indices correspondant a phi et theta
78 * respectivement.
79 * L'espace memoire correspondant a ce pointeur doit etre
80 * dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant
81 * l'appel a la routine.
82 *
83 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
84 * seul tableau, qui constitue une entree/sortie.
85 */
86
87/*
88 * $Id: cfrchebpii.C,v 1.4 2016/12/05 16:18:03 j_novak Exp $
89 * $Log: cfrchebpii.C,v $
90 * Revision 1.4 2016/12/05 16:18:03 j_novak
91 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
92 *
93 * Revision 1.3 2014/10/13 08:53:15 j_novak
94 * Lorene classes and functions now belong to the namespace Lorene.
95 *
96 * Revision 1.2 2014/10/06 15:18:44 j_novak
97 * Modified #include directives to use c++ syntax.
98 *
99 * Revision 1.1 2004/12/21 17:06:01 j_novak
100 * Added all files for using fftw3.
101 *
102 * Revision 1.1 2004/11/23 15:13:50 m_forot
103 * Added the bases for the cases without any equatorial symmetry
104 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
105 *
106 *
107 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cfrchebpii.C,v 1.4 2016/12/05 16:18:03 j_novak Exp $
108 *
109 */
110
111
112// headers du C
113#include <cassert>
114#include <cstdlib>
115
116// Prototypes of F77 subroutines
117#include "headcpp.h"
118#include "proto_f77.h"
119
120// Prototypage des sous-routines utilisees:
121namespace Lorene {
122int* facto_ini(int ) ;
123double* trigo_ini(int ) ;
124double* cheb_ini(const int) ;
125double* chebimp_ini(const int ) ;
126
127//*****************************************************************************
128
129void cfrchebpii(const int* deg, const int* dimf, double* ff, const int* dimc,
130 double* cf)
131
132{
133
134int i, j, k ;
135
136// Dimensions des tableaux ff et cf :
137 int n1f = dimf[0] ;
138 int n2f = dimf[1] ;
139 int n3f = dimf[2] ;
140 int n1c = dimc[0] ;
141 int n2c = dimc[1] ;
142 int n3c = dimc[2] ;
143
144// Nombres de degres de liberte en r :
145 int nr = deg[2] ;
146
147// Tests de dimension:
148 if (nr > n3f) {
149 cout << "cfrchebpii: nr > n3f : nr = " << nr << " , n3f = "
150 << n3f << endl ;
151 abort () ;
152 exit(-1) ;
153 }
154 if (nr > n3c) {
155 cout << "cfrchebpii: nr > n3c : nr = " << nr << " , n3c = "
156 << n3c << endl ;
157 abort () ;
158 exit(-1) ;
159 }
160 if (n1f > n1c) {
161 cout << "cfrchebpii: n1f > n1c : n1f = " << n1f << " , n1c = "
162 << n1c << endl ;
163 abort () ;
164 exit(-1) ;
165 }
166 if (n2f > n2c) {
167 cout << "cfrchebpii: n2f > n2c : n2f = " << n2f << " , n2c = "
168 << n2c << endl ;
169 abort () ;
170 exit(-1) ;
171 }
172
173// Nombre de points pour la FFT:
174 int nm1 = nr - 1;
175 int nm1s2 = nm1 / 2;
176
177// Recherche des tables pour la FFT:
178 int* facto = facto_ini(nm1) ;
179 double* trigo = trigo_ini(nm1) ;
180
181// Recherche de la table des sin(psi) :
182 double* sinp = cheb_ini(nr);
183
184// Recherche de la table des points de collocations x_k :
185 double* x = chebimp_ini(nr);
186
187 // tableau de travail G et t1
188 // (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
189 double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) );
190 double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
191
192// Parametres pour la routine FFT991
193 int jump = 1 ;
194 int inc = 1 ;
195 int lot = 1 ;
196 int isign = - 1 ;
197
198// boucle sur phi et theta
199
200 int n2n3f = n2f * n3f ;
201 int n2n3c = n2c * n3c ;
202
203/*
204 * Borne de la boucle sur phi:
205 * si n1f = 1, on effectue la boucle une fois seulement.
206 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
207 * j=n1f-1 et j=0 ne sont pas consideres car nuls).
208 */
209 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
210
211 for (j=0; j< borne_phi; j++) {
212
213 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
214
215 /************** Cas l impair ***************/
216
217 for (k=1; k<n2f; k+=2) {
218
219 int i0 = n2n3f * j + n3f * k ; // indice de depart
220 double* ff0 = ff + i0 ; // tableau des donnees a transformer
221
222 i0 = n2n3c * j + n3c * k ; // indice de depart
223 double* cf0 = cf + i0 ; // tableau resultat
224
225
226// Valeur en theta=0 de la partie antisymetrique de F, F_ :
227 double fmoins0 = 0.5 * ( ff0[nm1] - ff0[0] );
228
229// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
230//---------------------------------------------
231 for ( i = 1; i < nm1s2 ; i++ ) {
232// ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
233 int isym = nm1 - i ;
234// ... indice (dans le tableau ff0) du point x correspondant a theta
235 int ix = nm1 - i ;
236// ... indice (dans le tableau ff0) du point x correspondant a sym(theta)
237 int ixsym = nm1 - isym ;
238
239// ... F+(psi)
240 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
241// ... F_(psi) sin(psi)
242 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
243 g[i] = fp + fms ;
244 g[isym] = fp - fms ;
245 }
246//... cas particuliers:
247 g[0] = 0.5 * ( ff0[nm1] + ff0[0] );
248 g[nm1s2] = ff0[nm1s2];
249
250// Developpement de G en series de Fourier par une FFT
251//----------------------------------------------------
252
253 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
254
255// Coefficients pairs du developmt. de Tchebyshev de f
256//----------------------------------------------------
257// Ces coefficients sont egaux aux coefficients en cosinus du developpement
258// de G en series de Fourier (le facteur 2 vient de la normalisation
259// de fft991) :
260
261 cf0[0] = g[0] ;
262 for (i=2; i<nm1; i += 2 ) cf0[i] = 2. * g[i] ;
263 cf0[nm1] = g[nm1] ;
264
265// Coefficients impairs du developmt. de Tchebyshev de f
266//------------------------------------------------------
267// 1. Coef. c'_k (recurrence amorcee a partir de zero)
268// Le +4. en facteur de g[i] est du a la normalisation de fft991
269// (si fft991 donnait reellement les coef. en sinus, il faudrait le
270// remplacer par un -2.)
271 cf0[1] = 0 ;
272 double som = 0;
273 for ( i = 3; i < nr; i += 2 ) {
274 cf0[i] = cf0[i-2] + 4. * g[i] ;
275 som += cf0[i] ;
276 }
277
278// 2. Calcul de c_1 :
279 double c1 = ( fmoins0 - som ) / nm1s2 ;
280
281// 3. Coef. c_k avec k impair:
282 cf0[1] = c1 ;
283 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
284
285
286 } // fin de la boucle sur theta
287
288 /************** Cas l pair ***************/
289
290 for (k=0; k<n2f; k+=2) {
291 int i0 = n2n3f * j + n3f * k ; // indice de depart
292 double* ff0 = ff + i0 ; // tableau des donnees a transformer
293
294 i0 = n2n3c * j + n3c * k ; // indice de depart
295 double* cf0 = cf + i0 ; // tableau resultat
296
297// Multiplication de la fonction par x (pour la rendre paire)
298// NB: dans les commentaires qui suivent, on note h(x) = x f(x).
299// (Les valeurs de h dans l'espace des configurations sont stokees dans le
300// tableau cf0).
301 cf0[0] = 0 ;
302 for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
303
304// Valeur en psi=0 de la partie antisymetrique de F, F_ :
305 double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
306
307// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
308//---------------------------------------------
309 for ( i = 1; i < nm1s2 ; i++ ) {
310// ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
311 int isym = nm1 - i ;
312// ... indice (dans le tableau cf0) du point x correspondant a theta
313 int ix = nm1 - i ;
314// ... indice (dans le tableau cf0) du point x correspondant a sym(theta)
315 int ixsym = nm1 - isym ;
316
317// ... F+(psi)
318 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
319// ... F_(psi) sin(psi)
320 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
321 g[i] = fp + fms ;
322 g[isym] = fp - fms ;
323 }
324//... cas particuliers:
325 g[0] = 0.5 * ( cf0[nm1] + cf0[0] );
326 g[nm1s2] = cf0[nm1s2];
327
328// Developpement de G en series de Fourier par une FFT
329//----------------------------------------------------
330
331 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
332
333// Coefficients pairs du developmt. de Tchebyshev de h
334//----------------------------------------------------
335// Ces coefficients sont egaux aux coefficients en cosinus du developpement
336// de G en series de Fourier (le facteur 2. vient de la normalisation
337// de fft991; si fft991 donnait reellement les coef. en cosinus, il faudrait le
338// remplacer par un +1.) :
339
340 cf0[0] = g[0] ;
341 for (i=2; i<nm1; i += 2 ) cf0[i] = 2. * g[i] ;
342 cf0[nm1] = g[nm1] ;
343
344// Coefficients impairs du developmt. de Tchebyshev de h
345//------------------------------------------------------
346// 1. Coef. c'_k (recurrence amorcee a partir de zero)
347// Le +4. en facteur de g[i] est du a la normalisation de fft991
348// (si fft991 donnait reellement les coef. en sinus, il faudrait le
349// remplacer par un -2.)
350 cf0[1] = 0 ;
351 double som = 0;
352 for ( i = 3; i < nr; i += 2 ) {
353 cf0[i] = cf0[i-2] + 4. * g[i] ;
354 som += cf0[i] ;
355 }
356
357// 2. Calcul de c_1 :
358 double c1 = ( fmoins0 - som ) / nm1s2 ;
359
360// 3. Coef. c_k avec k impair:
361 cf0[1] = c1 ;
362 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
363
364// Coefficients de f en fonction de ceux de h
365//-------------------------------------------
366
367 cf0[0] = 2* cf0[0] ;
368 for (i=1; i<nm1; i++) {
369 cf0[i] = 2 * cf0[i] - cf0[i-1] ;
370 }
371 cf0[nm1] = 0 ;
372
373
374 } // fin de la boucle sur theta
375 } // fin de la boucle sur phi
376
377 // Menage
378 free (t1) ;
379 free (g) ;
380
381}
382}
383
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