LORENE
FFTW3/citcosi.C
1/*
2 * Copyright (c) 1999-2001 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 en cos((2*l+1)*theta) inverse sur le deuxieme indice (theta)
27 * d'un tableau 3-D representant une fonction antisymetrique par rapport
28 * au plan z=0.
29 * Utilise la bibliotheque fftw.
30 *
31 * Entree:
32 * -------
33 * int* deg : tableau du nombre effectif de degres de liberte dans chacune
34 * des 3 dimensions: le nombre de points de collocation
35 * en theta est nt = deg[1] et doit etre de la forme
36 * nt = 2*p + 1
37 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
38 * dimensions.
39 * On doit avoir dimc[1] >= deg[1] = nt.
40 * NB: pour dimc[0] = 1 (un seul point en phi), la transformation
41 * est bien effectuee.
42 * pour dimc[0] > 1 (plus d'un point en phi), la
43 * transformation n'est effectuee que pour les indices (en phi)
44 * j != 1 et j != dimc[0]-1 (cf. commentaires sur borne_phi).
45 *
46 * double* cf : tableau des coefficients c_l de la fonction definis
47 * comme suit (a r et phi fixes)
48 *
49 * f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) .
50 *
51 * L'espace memoire correspondant a ce
52 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
53 * avoir ete alloue avant l'appel a la routine.
54 * Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
55 * le tableau cf comme suit
56 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
57 * ou j et k sont les indices correspondant a
58 * phi et r respectivement.
59 *
60 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
61 * dimensions.
62 * On doit avoir dimf[1] >= deg[1] = nt.
63 *
64 * Sortie:
65 * -------
66 * double* ff : tableau des valeurs de la fonction aux nt points de
67 * de collocation
68 *
69 * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
70 *
71 * L'espace memoire correspondant a ce
72 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
73 * avoir ete alloue avant l'appel a la routine.
74 * Les valeurs de la fonction sont stokees
75 * dans le tableau ff comme suit
76 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
77 * ou j et k sont les indices correspondant a
78 * phi et r respectivement.
79 *
80 * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
81 * seul tableau, qui constitue une entree/sortie.
82 *
83 */
84
85/*
86 * $Id: citcosi.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
87 * $Log: citcosi.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:20 j_novak
92 * Lorene classes and functions now belong to the namespace Lorene.
93 *
94 * Revision 1.2 2014/10/06 15:18:49 j_novak
95 * Modified #include directives to use c++ syntax.
96 *
97 * Revision 1.1 2004/12/21 17:06:03 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:53 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:40 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:42:54 hyc
118 * *** empty log message ***
119 *
120 *
121 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcosi.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
122 *
123 */
124
125
126// headers du C
127#include <cstdlib>
128#include <fftw3.h>
129
130//Lorene prototypes
131#include "tbl.h"
132
133// Prototypage des sous-routines utilisees:
134namespace Lorene {
135fftw_plan back_fft(int, Tbl*&) ;
136double* cheb_ini(const int) ;
137double* chebimp_ini(const int ) ;
138//*****************************************************************************
139
140void citcosi(const int* deg, const int* dimc, double* cf, const int* dimf,
141 double* ff)
142{
143
144int i, j, k ;
145
146// Dimensions des tableaux ff et cf :
147 int n1f = dimf[0] ;
148 int n2f = dimf[1] ;
149 int n3f = dimf[2] ;
150 int n1c = dimc[0] ;
151 int n2c = dimc[1] ;
152 int n3c = dimc[2] ;
153
154// Nombres de degres de liberte en theta :
155 int nt = deg[1] ;
156
157// Tests de dimension:
158 if (nt > n2f) {
159 cout << "citcosi: nt > n2f : nt = " << nt << " , n2f = "
160 << n2f << endl ;
161 abort () ;
162 }
163 if (nt > n2c) {
164 cout << "citcosi: nt > n2c : nt = " << nt << " , n2c = "
165 << n2c << endl ;
166 abort () ;
167 }
168 if ( (n1f > 1) && (n1c > n1f) ) {
169 cout << "citcosi: n1c > n1f : n1c = " << n1c << " , n1f = "
170 << n1f << endl ;
171 abort () ;
172 }
173 if (n3c > n3f) {
174 cout << "citcosi: n3c > n3f : n3c = " << n3c << " , n3f = "
175 << n3f << endl ;
176 abort () ;
177 }
178
179// Nombre de points pour la FFT:
180 int nm1 = nt - 1;
181 int nm1s2 = nm1 / 2;
182
183// Recherche des tables pour la FFT:
184 Tbl* pg = 0x0 ;
185 fftw_plan p = back_fft(nm1, pg) ;
186 Tbl& g = *pg ;
187 double* t1 = new double[nt] ;
188
189// Recherche de la table des sin(psi) :
190 double* sinp = cheb_ini(nt) ;
191
192// Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}):
193 double* x = chebimp_ini(nt) ;
194
195// boucle sur phi et r
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 = n1c-2 en sautant j = 1 (les coefficients
204 * j=n1c-1 et j=0 ne sont pas consideres car nuls).
205 */
206 int borne_phi = n1c-1 ;
207 if (n1f == 1) borne_phi = 1 ;
208
209 for (j=0; j< borne_phi; j++) {
210
211 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
212
213 for (k=0; k<n3c; k++) {
214
215 int i0 = n2n3c * j + k ; // indice de depart
216 double* cf0 = cf + i0 ; // tableau des donnees a transformer
217
218 i0 = n2n3f * j + k ; // indice de depart
219 double* ff0 = ff + i0 ; // tableau resultat
220
221// Calcul des coefficients du developpement en T_{2i}(x) de la fonction
222// h(x) := x f(x) (x=cos(theta)) a partir des coefficients de f
223// (resultat stoke dans le tableau t1 :
224 t1[0] = .5 * cf0[0] ;
225 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
226 t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
227
228/*
229 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
230 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
231 */
232
233// Calcul des coefficients de Fourier de la fonction
234// G(psi) = F+(psi) + F_(psi) sin(psi)
235// en fonction des coefficients en cos(2l theta) de f:
236
237// Coefficients impairs de G
238//--------------------------
239
240 double c1 = t1[1] ;
241
242 double som = 0;
243 ff0[n3f] = 0 ;
244 for ( i = 3; i < nt; i += 2 ) {
245 ff0[ n3f*i ] = t1[i] - c1 ;
246 som += ff0[ n3f*i ] ;
247 }
248
249// Valeur en psi=0 de la partie antisymetrique de F, F_ :
250 double fmoins0 = nm1s2 * c1 + som ;
251
252// Coef. impairs de G
253// NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
254// donnait exactement les coef. des sinus, ce facteur serait -0.5.
255 for ( i = 3; i < nt; i += 2 ) {
256 g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
257 }
258
259
260// Coefficients pairs de G
261//------------------------
262// Ces coefficients sont egaux aux coefficients pairs du developpement de
263// f.
264// NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
265// donnait exactement les coef. des cosinus, ce facteur serait 1.
266
267 g.set(0) = t1[0] ;
268 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
269 g.set(nm1s2) = t1[nm1] ;
270
271// Transformation de Fourier inverse de G
272//---------------------------------------
273
274// FFT inverse
275 fftw_execute(p) ;
276
277// Valeurs de f deduites de celles de G
278//-------------------------------------
279
280 for ( i = 1; i < nm1s2 ; i++ ) {
281// ... indice du pt symetrique de psi par rapport a pi/2:
282 int isym = nm1 - i ;
283
284 double fp = 0.5 * ( g(i) + g(isym) ) ;
285 double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
286 ff0[ n3f*i ] = ( fp + fm ) / x[isym] ;
287 ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
288 }
289
290//... cas particuliers:
291 ff0[0] = g(0) + fmoins0 ;
292 ff0[ n3f*nm1 ] = 0 ;
293 ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
294
295
296 } // fin de la boucle sur r
297 } // fin de la boucle sur phi
298
299 delete [] t1 ;
300
301}
302}
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