LORENE
FFTW3/cfrchebpip.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 bibliotheque fftw.
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 + 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 pair f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) ,
71 * Si l impair 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: cfrchebpip.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
89 * $Log: cfrchebpip.C,v $
90 * Revision 1.4 2016/12/05 16:18:04 j_novak
91 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
92 *
93 * Revision 1.3 2014/10/13 08:53:18 j_novak
94 * Lorene classes and functions now belong to the namespace Lorene.
95 *
96 * Revision 1.2 2014/10/06 15:18:48 j_novak
97 * Modified #include directives to use c++ syntax.
98 *
99 * Revision 1.1 2004/12/21 17:06:02 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/FFTW3/cfrchebpip.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
108 *
109 */
110
111
112// headers du C
113#include <cstdlib>
114#include <fftw3.h>
115
116//Lorene prototypes
117#include "tbl.h"
118
119// Prototypage des sous-routines utilisees:
120namespace Lorene {
121fftw_plan prepare_fft(int, Tbl*&) ;
122double* cheb_ini(const int) ;
123double* chebimp_ini(const int ) ;
124
125//*****************************************************************************
126
127void cfrchebpip(const int* deg, const int* dimf, double* ff, const int* dimc,
128 double* cf)
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 > n3f) {
147 cout << "cfrchebpip: nr > n3f : nr = " << nr << " , n3f = "
148 << n3f << endl ;
149 abort () ;
150 exit(-1) ;
151 }
152 if (nr > n3c) {
153 cout << "cfrchebpip: nr > n3c : nr = " << nr << " , n3c = "
154 << n3c << endl ;
155 abort () ;
156 exit(-1) ;
157 }
158 if (n1f > n1c) {
159 cout << "cfrchebpip: n1f > n1c : n1f = " << n1f << " , n1c = "
160 << n1c << endl ;
161 abort () ;
162 exit(-1) ;
163 }
164 if (n2f > n2c) {
165 cout << "cfrchebpip: n2f > n2c : n2f = " << n2f << " , n2c = "
166 << n2c << 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 Tbl* pg = 0x0 ;
177 fftw_plan p = prepare_fft(nm1, pg) ;
178 Tbl& g = *pg ;
179
180// Recherche de la table des sin(psi) :
181 double* sinp = cheb_ini(nr);
182
183// Recherche de la table des points de collocations x_k :
184 double* x = chebimp_ini(nr);
185
186// boucle sur phi et theta
187
188 int n2n3f = n2f * n3f ;
189 int n2n3c = n2c * n3c ;
190
191/*
192 * Borne de la boucle sur phi:
193 * si n1f = 1, on effectue la boucle une fois seulement.
194 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
195 * j=n1f-1 et j=0 ne sont pas consideres car nuls).
196 */
197 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
198
199 for (j=0; j< borne_phi; j++) {
200
201 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
202
203 /************** Cas l pair ***************/
204
205 for (k=0; k<n2f; k+=2) {
206
207 int i0 = n2n3f * j + n3f * k ; // indice de depart
208 double* ff0 = ff + i0 ; // tableau des donnees a transformer
209
210 i0 = n2n3c * j + n3c * k ; // indice de depart
211 double* cf0 = cf + i0 ; // tableau resultat
212
213
214// Valeur en theta=0 de la partie antisymetrique de F, F_ :
215 double fmoins0 = 0.5 * ( ff0[nm1] - ff0[0] );
216
217// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
218//---------------------------------------------
219 for ( i = 1; i < nm1s2 ; i++ ) {
220// ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
221 int isym = nm1 - i ;
222// ... indice (dans le tableau ff0) du point x correspondant a theta
223 int ix = nm1 - i ;
224// ... indice (dans le tableau ff0) du point x correspondant a sym(theta)
225 int ixsym = nm1 - isym ;
226
227// ... F+(psi)
228 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
229// ... F_(psi) sin(psi)
230 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
231 g.set(i) = fp + fms ;
232 g.set(isym) = fp - fms ;
233 }
234//... cas particuliers:
235 g.set(0) = 0.5 * ( ff0[nm1] + ff0[0] );
236 g.set(nm1s2) = ff0[nm1s2];
237
238// Developpement de G en series de Fourier par une FFT
239//----------------------------------------------------
240
241 fftw_execute(p) ;
242
243// Coefficients pairs du developmt. de Tchebyshev de f
244//----------------------------------------------------
245// Ces coefficients sont egaux aux coefficients en cosinus du developpement
246// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
247// de fftw) :
248
249 double fac = 2./double(nm1) ;
250 cf0[0] = g(0) / double(nm1) ;
251 for (i=2; i<nm1; i += 2) cf0[i] = fac*g(i/2) ;
252 cf0[nm1] = g(nm1s2) / double(nm1) ;
253
254// Coefficients impairs du developmt. de Tchebyshev de f
255//------------------------------------------------------
256// 1. Coef. c'_k (recurrence amorcee a partir de zero)
257// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
258// (si fftw donnait reellement les coef. en sinus, il faudrait le
259// remplacer par un -2.)
260 fac *= 2. ;
261 cf0[1] = 0. ;
262 double som = 0.;
263 for ( i = 3; i < nr; i += 2 ) {
264 cf0[i] = cf0[i-2] + fac * g(nm1-i/2) ;
265 som += cf0[i] ;
266 }
267
268// 2. Calcul de c_1 :
269 double c1 = ( fmoins0 - som ) / nm1s2 ;
270
271// 3. Coef. c_k avec k impair:
272 cf0[1] = c1 ;
273 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
274
275
276 } // fin de la boucle sur theta
277
278 /************** Cas l impair ***************/
279
280 for (k=1; k<n2f; k+=2) {
281 int i0 = n2n3f * j + n3f * k ; // indice de depart
282 double* ff0 = ff + i0 ; // tableau des donnees a transformer
283
284 i0 = n2n3c * j + n3c * k ; // indice de depart
285 double* cf0 = cf + i0 ; // tableau resultat
286
287// Multiplication de la fonction par x (pour la rendre paire)
288// NB: dans les commentaires qui suivent, on note h(x) = x f(x).
289// (Les valeurs de h dans l'espace des configurations sont stokees dans le
290// tableau cf0).
291 cf0[0] = 0 ;
292 for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
293
294
295// Valeur en psi=0 de la partie antisymetrique de F, F_ :
296 double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
297
298// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
299//---------------------------------------------
300 for ( i = 1; i < nm1s2 ; i++ ) {
301// ... indice (dans le tableau g) du pt symetrique de theta par rapport a pi/2:
302 int isym = nm1 - i ;
303// ... indice (dans le tableau cf0) du point x correspondant a theta
304 int ix = nm1 - i ;
305// ... indice (dans le tableau cf0) du point x correspondant a sym(theta)
306 int ixsym = nm1 - isym ;
307
308// ... F+(psi)
309 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
310// ... F_(psi) sin(psi)
311 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
312 g.set(i) = fp + fms ;
313 g.set(isym) = fp - fms ;
314 }
315//... cas particuliers:
316 g.set(0) = 0.5 * ( cf0[nm1] + cf0[0] );
317 g.set(nm1s2) = cf0[nm1s2];
318
319// Developpement de G en series de Fourier par une FFT
320//----------------------------------------------------
321
322 fftw_execute(p) ;
323
324// Coefficients pairs du developmt. de Tchebyshev de h
325//----------------------------------------------------
326// Ces coefficients sont egaux aux coefficients en cosinus du developpement
327// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
328// de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
329// remplacer par un +1.) :
330
331 double fac = 2./double(nm1) ;
332 cf0[0] = g(0) / double(nm1) ;
333 for (i=2; i<nm1; i += 2 ) cf0[i] = fac * g(i/2) ;
334 cf0[nm1] = g(nm1s2) / double(nm1) ;
335
336// Coefficients impairs du developmt. de Tchebyshev de h
337//------------------------------------------------------
338// 1. Coef. c'_k (recurrence amorcee a partir de zero)
339// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
340// (si fftw donnait reellement les coef. en sinus, il faudrait le
341// remplacer par un -2.)
342 fac *= 2 ;
343 cf0[1] = 0 ;
344 double som = 0;
345 for ( i = 3; i < nr; i += 2 ) {
346 cf0[i] = cf0[i-2] + fac * g(nm1 - i/2) ;
347 som += cf0[i] ;
348 }
349
350// 2. Calcul de c_1 :
351 double c1 = ( fmoins0 - som ) / nm1s2 ;
352
353// 3. Coef. c_k avec k impair:
354 cf0[1] = c1 ;
355 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
356
357// Coefficients de f en fonction de ceux de h
358//-------------------------------------------
359
360 cf0[0] = 2* cf0[0] ;
361 for (i=1; i<nm1; i++) {
362 cf0[i] = 2 * cf0[i] - cf0[i-1] ;
363 }
364 cf0[nm1] = 0 ;
365
366
367 } // fin de la boucle sur theta
368 } // fin de la boucle sur phi
369
370}
371}
Basic array class.
Definition tbl.h:161
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