LORENE
FFTW3/cftcosp.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 en cos(2*l*theta) sur le deuxieme indice (theta)
28 * d'un tableau 3-D representant une fonction symetrique par rapport
29 * au plan z=0.
30 * Utilise la bibliotheque fftw.
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 theta est nt = deg[1] et doit etre de la forme
37 * nt = 2*p + 1
38 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
39 * dimensions.
40 * On doit avoir dimf[1] >= deg[1] = nt.
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 nt points de
48 * de collocation
49 *
50 * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
51 *
52 * L'espace memoire correspondant a ce
53 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
54 * etre alloue avant l'appel a la routine.
55 * Les valeurs de la fonction doivent etre stokees
56 * dans le tableau ff comme suit
57 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
58 * ou j et k sont les indices correspondant a
59 * phi et r respectivement.
60 *
61 * int* dimc : tableau du nombre d'elements de cc dans chacune des trois
62 * dimensions.
63 * On doit avoir dimc[1] >= deg[1] = nt.
64 * Sortie:
65 * -------
66 * double* cf : tableau des coefficients c_l de la fonction definis
67 * comme suit (a r et phi fixes)
68 *
69 * f(theta) = som_{l=0}^{nt-1} c_l cos( 2 l theta ) .
70 *
71 * L'espace memoire correspondant a ce
72 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
73 * etre alloue avant l'appel a la routine.
74 * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
75 * le tableau cf comme suit
76 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
77 * ou j et k sont les indices correspondant a
78 * phi et r respectivement.
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/*
86 * $Id: cftcosp.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
87 * $Log: cftcosp.C,v $
88 * Revision 1.4 2016/12/05 16:18:05 j_novak
89 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
90 *
91 * Revision 1.3 2014/10/13 08:53:19 j_novak
92 * Lorene classes and functions now belong to the namespace Lorene.
93 *
94 * Revision 1.2 2014/10/06 15:18:48 j_novak
95 * Modified #include directives to use c++ syntax.
96 *
97 * Revision 1.1 2004/12/21 17:06:02 j_novak
98 * Added all files for using fftw3.
99 *
100 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
101 * Suppressed the directive #include <malloc.h> for malloc is defined
102 * in <stdlib.h>
103 *
104 * Revision 1.3 2002/10/16 14:36:44 j_novak
105 * Reorganization of #include instructions of standard C++, in order to
106 * use experimental version 3 of gcc.
107 *
108 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
109 * Modification of declaration of Fortran 77 prototypes for
110 * a better portability (in particular on IBM AIX systems):
111 * All Fortran subroutine names are now written F77_* and are
112 * defined in the new file C++/Include/proto_f77.h.
113 *
114 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
115 * LORENE
116 *
117 * Revision 2.0 1999/02/22 15:47:50 hyc
118 * *** empty log message ***
119 *
120 *
121 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcosp.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
122 *
123 */
124
125// headers du C
126#include <cstdlib>
127#include <fftw3.h>
128
129//Lorene prototypes
130#include "tbl.h"
131
132// Prototypage des sous-routines utilisees:
133namespace Lorene {
134fftw_plan prepare_fft(int, Tbl*&) ;
135double* cheb_ini(const int) ;
136
137//*****************************************************************************
138
139void cftcosp(const int* deg, const int* dimf, double* ff, const int* dimc,
140 double* cf)
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// Nombre de degres de liberte en theta :
154 int nt = deg[1] ;
155
156// Tests de dimension:
157 if (nt > n2f) {
158 cout << "cftcosp: nt > n2f : nt = " << nt << " , n2f = "
159 << n2f << endl ;
160 abort () ;
161 exit(-1) ;
162 }
163 if (nt > n2c) {
164 cout << "cftcosp: nt > n2c : nt = " << nt << " , n2c = "
165 << n2c << endl ;
166 abort () ;
167 exit(-1) ;
168 }
169 if (n1f > n1c) {
170 cout << "cftcosp: n1f > n1c : n1f = " << n1f << " , n1c = "
171 << n1c << endl ;
172 abort () ;
173 exit(-1) ;
174 }
175 if (n3f > n3c) {
176 cout << "cftcosp: n3f > n3c : n3f = " << n3f << " , n3c = "
177 << n3c << endl ;
178 abort () ;
179 exit(-1) ;
180 }
181
182// Nombre de points pour la FFT:
183 int nm1 = nt - 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(nt);
193
194// boucle sur phi et r (les boucles vont resp. de 0 a max(dimf[0]-2,0) et
195// 0 a dimf[2]-1 )
196
197 int n2n3f = n2f * n3f ;
198 int n2n3c = n2c * n3c ;
199
200/*
201 * Borne de la boucle sur phi:
202 * si n1f = 1, on effectue la boucle une fois seulement.
203 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
204 * j=n1f-1 et j=0 ne sont pas consideres car nuls).
205 */
206 int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
207
208 for (j=0; j< borne_phi; j++) {
209
210 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
211
212 for (k=0; k<n3f; k++) {
213
214 int i0 = n2n3f * j + k ; // indice de depart
215 double* ff0 = ff + i0 ; // tableau des donnees a transformer
216
217 i0 = n2n3c * j + k ; // indice de depart
218 double* cf0 = cf + i0 ; // tableau resultat
219
220/*
221 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
222 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
223 */
224
225// Valeur en psi=0 de la partie antisymetrique de F, F_ :
226 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
227
228// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
229//---------------------------------------------
230 for ( i = 1; i < nm1s2 ; i++ ) {
231// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
232 int isym = nm1 - i ;
233// ... indice (dans le tableau ff0) du point theta correspondant a psi
234 int ix = n3f * i ;
235// ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
236 int ixsym = n3f * isym ;
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[0] + ff0[ n3f*nm1 ] );
246 g.set(nm1s2) = ff0[ n3f*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. cos(2l theta) 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[n3c*i] = fac*g(i/2) ;
262 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
263
264// Coefficients impairs du developmt. en cos(2l theta) 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[n3c] = 0 ;
272 double som = 0;
273 for ( i = 3; i < nt; i += 2 ) {
274 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
275 som += cf0[n3c*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[n3c] = c1 ;
283 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
284
285
286 } // fin de la boucle sur r
287 } // fin de la boucle sur phi
288
289}
290}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord sinp
Definition map.h:735