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