LORENE
FFTW3/cfrcheb.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 * Transformation de Tchebyshev (cas fin) sur le troisieme indice (indice
27 * correspondant a r) d'un tableau 3-D
28 * par le biais de la bibliotheque fftw.
29 *
30 * Entree:
31 * -------
32 * int* deg : tableau du nombre effectif de degres de liberte dans chacune
33 * des 3 dimensions: le nombre de points de collocation
34 * en r est nr = deg[2] et doit etre de la forme
35 * nr = 2*p + 1
36 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
37 * dimensions.
38 * On doit avoir dimf[2] >= deg[2] = nr.
39 * NB: pour dimf[0] = 1 (un seul point en phi), la transformation
40 * est bien effectuee.
41 * pour dimf[0] > 1 (plus d'un point en phi), la
42 * transformation n'est effectuee que pour les indices (en phi)
43 * j != 1 et j != dimf[0]-1 (cf. commentaires sur borne_phi).
44 *
45 * double* ff : tableau des valeurs de la fonction aux nr points de
46 * de collocation
47 *
48 * x_i = - cos( pi i/(nr-1) ) 0 <= i <= nr-1
49 *
50 * Les valeurs de la fonction doivent etre stokees dans le
51 * tableau ff comme suit
52 * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ]
53 * ou j et k sont les indices correspondant a phi et theta
54 * respectivement.
55 * L'espace memoire correspondant a ce pointeur doit etre
56 * dimf[0]*dimf[1]*dimf[2] et doit avoir ete alloue avant
57 * l'appel a la routine.
58 *
59 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
60 * dimensions.
61 * On doit avoir dimc[2] >= deg[2] = nr.
62 *
63 * Sortie:
64 * -------
65 * double* cf : tableau des coefficients c_i de la fonction definis
66 * comme suit (a theta et phi fixes)
67 *
68 * f(x) = som_{i=0}^{nr-1} c_i T_i(x) ,
69 *
70 * ou T_i(x) designe le polynome de Tchebyshev de degre i.
71 * Les coefficients c_i (0 <= i <= nr-1) sont stokes dans
72 * le tableau cf comme suit
73 * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ]
74 * ou j et k sont les indices correspondant a phi et theta
75 * respectivement.
76 * L'espace memoire correspondant au pointeur cf doit etre
77 * dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant
78 * l'appel a la routine.
79 *
80 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
81 * seul tableau, qui constitue une entree/sortie.
82 */
83
84/*
85 * $Id: cfrcheb.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
86 * $Log: cfrcheb.C,v $
87 * Revision 1.4 2016/12/05 16:18:04 j_novak
88 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
89 *
90 * Revision 1.3 2014/10/13 08:53:18 j_novak
91 * Lorene classes and functions now belong to the namespace Lorene.
92 *
93 * Revision 1.2 2014/10/06 15:18:47 j_novak
94 * Modified #include directives to use c++ syntax.
95 *
96 * Revision 1.1 2004/12/21 17:06:02 j_novak
97 * Added all files for using fftw3.
98 *
99 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
100 * Suppressed the directive #include <malloc.h> for malloc is defined
101 * in <stdlib.h>
102 *
103 * Revision 1.3 2002/10/16 14:36:43 j_novak
104 * Reorganization of #include instructions of standard C++, in order to
105 * use experimental version 3 of gcc.
106 *
107 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
108 * Modification of declaration of Fortran 77 prototypes for
109 * a better portability (in particular on IBM AIX systems):
110 * All Fortran subroutine names are now written F77_* and are
111 * defined in the new file C++/Include/proto_f77.h.
112 *
113 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
114 * LORENE
115 *
116 * Revision 2.0 1999/02/22 15:48:47 hyc
117 * *** empty log message ***
118 *
119 *
120 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfrcheb.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
121 *
122 */
123
124
125// headers du C
126#include <cstdlib>
127#include <cassert>
128#include <fftw3.h>
129
130//Lorene prototypes
131#include "tbl.h"
132
133// Prototypage des sous-routines utilisees:
134namespace Lorene {
135fftw_plan prepare_fft(int, Tbl*&) ;
136double* cheb_ini(const int) ;
137
138//*****************************************************************************
139
140void cfrcheb(const int* deg, const int* dimf, double* ff, const int* dimc,
141 double* cf)
142
143{
144
145int i, j, k ;
146
147// Dimensions des tableaux ff et cf :
148 int n1f = dimf[0] ;
149 int n2f = dimf[1] ;
150 int n3f = dimf[2] ;
151 int n1c = dimc[0] ;
152 int n2c = dimc[1] ;
153 int n3c = dimc[2] ;
154
155// Nombres de degres de liberte en r :
156 int nr = deg[2] ;
157
158// Tests de dimension:
159 if (nr > n3f) {
160 cout << "cfrcheb: nr > n3f : nr = " << nr << " , n3f = "
161 << n3f << endl ;
162 abort () ;
163 exit(-1) ;
164 }
165 if (nr > n3c) {
166 cout << "cfrcheb: nr > n3c : nr = " << nr << " , n3c = "
167 << n3c << endl ;
168 abort () ;
169 exit(-1) ;
170 }
171 if (n1f > n1c) {
172 cout << "cfrcheb: n1f > n1c : n1f = " << n1f << " , n1c = "
173 << n1c << endl ;
174 abort () ;
175 exit(-1) ;
176 }
177 if (n2f > n2c) {
178 cout << "cfrcheb: n2f > n2c : n2f = " << n2f << " , n2c = "
179 << n2c << endl ;
180 abort () ;
181 exit(-1) ;
182 }
183
184
185// Nombre de points pour la FFT:
186 int nm1 = nr - 1;
187 int nm1s2 = nm1 / 2;
188
189// Recherche des tables pour la FFT:
190 Tbl* pg = 0x0 ;
191 fftw_plan p = prepare_fft(nm1, pg) ;
192 Tbl& g = *pg ;
193
194// Recherche de la table des sin(psi) :
195 double* sinp = cheb_ini(nr);
196
197// boucle sur phi et theta
198
199 int n2n3f = n2f * n3f ;
200 int n2n3c = n2c * n3c ;
201
202/*
203 * Borne de la boucle sur phi:
204 * si n1f = 1, on effectue la boucle une fois seulement.
205 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
206 * j=n1f-1 et j=0 ne sont pas consideres car nuls).
207 */
208 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
209
210 for (j=0; j< borne_phi; j++) {
211
212 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
213
214 for (k=0; k<n2f; k++) {
215
216 int i0 = n2n3f * j + n3f * k ; // indice de depart
217 double* ff0 = ff + i0 ; // tableau des donnees a transformer
218
219 i0 = n2n3c * j + n3c * k ; // indice de depart
220 double* cf0 = cf + i0 ; // tableau resultat
221
222/*
223 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
224 * reliee a x par x = - cos(psi) et F(psi) = f(x(psi)).
225 */
226
227// Valeur en psi=0 de la partie antisymetrique de F, F_ :
228 double fmoins0 = 0.5 * ( ff0[0] - ff0[nm1] );
229
230// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
231//---------------------------------------------
232 for ( i = 1; i < nm1s2 ; i++ ) {
233// ... indice du pt symetrique de psi par rapport a pi/2:
234 int isym = nm1 - i ;
235// ... F+(psi)
236 double fp = 0.5 * ( ff0[i] + ff0[isym] ) ;
237// ... F_(psi) sin(psi)
238 double fms = 0.5 * ( ff0[i] - ff0[isym] ) * sinp[i] ;
239 g.set(i) = fp + fms ;
240 g.set(isym) = fp - fms ;
241 }
242//... cas particuliers:
243 g.set(0) = 0.5 * ( ff0[0] + ff0[nm1] );
244 g.set(nm1s2) = ff0[nm1s2];
245
246// Developpement de G en series de Fourier par une FFT
247//----------------------------------------------------
248
249 fftw_execute(p) ;
250
251// Coefficients pairs du developmt. de Tchebyshev de f
252//----------------------------------------------------
253// Ces coefficients sont egaux aux coefficients en cosinus du developpement
254// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
255// de fftw) :
256
257 double fac = 2./double(nm1) ;
258 cf0[0] = g(0) / double(nm1) ;
259 for (i=2; i<nm1; i += 2) cf0[i] = fac*g(i/2) ;
260 cf0[nm1] = g(nm1s2) / double(nm1) ;
261
262// Coefficients impairs du developmt. de Tchebyshev de f
263//------------------------------------------------------
264// 1. Coef. c'_k (recurrence amorcee a partir de zero):
265// NB: Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
266// (si fftw donnait reellement les coef. en sinus, il faudrait le
267// remplacer par un +2.)
268 fac *= -2. ;
269 cf0[1] = 0 ;
270 double som = 0;
271 for ( i = 3; i < nr; i += 2 ) {
272 cf0[i] = cf0[i-2] + fac * g(nm1-i/2) ;
273 som += cf0[i] ;
274 }
275
276// 2. Calcul de c_1 :
277 double c1 = - ( fmoins0 + som ) / nm1s2 ;
278
279// 3. Coef. c_k avec k impair:
280 cf0[1] = c1 ;
281 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
282
283
284 } // fin de la boucle sur theta
285 } // fin de la boucle sur phi
286
287}
288}
289
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord sinp
Definition map.h:735