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
23char citcosi_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcosi.C,v 1.3 2014/10/13 08:53:20 j_novak Exp $" ;
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.3 2014/10/13 08:53:20 j_novak Exp $
87 * $Log: citcosi.C,v $
88 * Revision 1.3 2014/10/13 08:53:20 j_novak
89 * Lorene classes and functions now belong to the namespace Lorene.
90 *
91 * Revision 1.2 2014/10/06 15:18:49 j_novak
92 * Modified #include directives to use c++ syntax.
93 *
94 * Revision 1.1 2004/12/21 17:06:03 j_novak
95 * Added all files for using fftw3.
96 *
97 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
98 * Suppressed the directive #include <malloc.h> for malloc is defined
99 * in <stdlib.h>
100 *
101 * Revision 1.3 2002/10/16 14:36:53 j_novak
102 * Reorganization of #include instructions of standard C++, in order to
103 * use experimental version 3 of gcc.
104 *
105 * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
106 * Modification of declaration of Fortran 77 prototypes for
107 * a better portability (in particular on IBM AIX systems):
108 * All Fortran subroutine names are now written F77_* and are
109 * defined in the new file C++/Include/proto_f77.h.
110 *
111 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
112 * LORENE
113 *
114 * Revision 2.0 1999/02/22 15:42:54 hyc
115 * *** empty log message ***
116 *
117 *
118 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcosi.C,v 1.3 2014/10/13 08:53:20 j_novak Exp $
119 *
120 */
121
122
123// headers du C
124#include <cstdlib>
125#include <fftw3.h>
126
127//Lorene prototypes
128#include "tbl.h"
129
130// Prototypage des sous-routines utilisees:
131namespace Lorene {
132fftw_plan back_fft(int, Tbl*&) ;
133double* cheb_ini(const int) ;
134double* chebimp_ini(const int ) ;
135//*****************************************************************************
136
137void citcosi(const int* deg, const int* dimc, double* cf, const int* dimf,
138 double* ff)
139{
140
141int i, j, k ;
142
143// Dimensions des tableaux ff et cf :
144 int n1f = dimf[0] ;
145 int n2f = dimf[1] ;
146 int n3f = dimf[2] ;
147 int n1c = dimc[0] ;
148 int n2c = dimc[1] ;
149 int n3c = dimc[2] ;
150
151// Nombres de degres de liberte en theta :
152 int nt = deg[1] ;
153
154// Tests de dimension:
155 if (nt > n2f) {
156 cout << "citcosi: nt > n2f : nt = " << nt << " , n2f = "
157 << n2f << endl ;
158 abort () ;
159 }
160 if (nt > n2c) {
161 cout << "citcosi: nt > n2c : nt = " << nt << " , n2c = "
162 << n2c << endl ;
163 abort () ;
164 }
165 if ( (n1f > 1) && (n1c > n1f) ) {
166 cout << "citcosi: n1c > n1f : n1c = " << n1c << " , n1f = "
167 << n1f << endl ;
168 abort () ;
169 }
170 if (n3c > n3f) {
171 cout << "citcosi: n3c > n3f : n3c = " << n3c << " , n3f = "
172 << n3f << endl ;
173 abort () ;
174 }
175
176// Nombre de points pour la FFT:
177 int nm1 = nt - 1;
178 int nm1s2 = nm1 / 2;
179
180// Recherche des tables pour la FFT:
181 Tbl* pg = 0x0 ;
182 fftw_plan p = back_fft(nm1, pg) ;
183 Tbl& g = *pg ;
184 double* t1 = new double[nt] ;
185
186// Recherche de la table des sin(psi) :
187 double* sinp = cheb_ini(nt) ;
188
189// Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}):
190 double* x = chebimp_ini(nt) ;
191
192// boucle sur phi et r
193
194 int n2n3f = n2f * n3f ;
195 int n2n3c = n2c * n3c ;
196
197/*
198 * Borne de la boucle sur phi:
199 * si n1f = 1, on effectue la boucle une fois seulement.
200 * si n1f > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
201 * j=n1c-1 et j=0 ne sont pas consideres car nuls).
202 */
203 int borne_phi = n1c-1 ;
204 if (n1f == 1) borne_phi = 1 ;
205
206 for (j=0; j< borne_phi; j++) {
207
208 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
209
210 for (k=0; k<n3c; k++) {
211
212 int i0 = n2n3c * j + k ; // indice de depart
213 double* cf0 = cf + i0 ; // tableau des donnees a transformer
214
215 i0 = n2n3f * j + k ; // indice de depart
216 double* ff0 = ff + i0 ; // tableau resultat
217
218// Calcul des coefficients du developpement en T_{2i}(x) de la fonction
219// h(x) := x f(x) (x=cos(theta)) a partir des coefficients de f
220// (resultat stoke dans le tableau t1 :
221 t1[0] = .5 * cf0[0] ;
222 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
223 t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
224
225/*
226 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
227 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
228 */
229
230// Calcul des coefficients de Fourier de la fonction
231// G(psi) = F+(psi) + F_(psi) sin(psi)
232// en fonction des coefficients en cos(2l theta) de f:
233
234// Coefficients impairs de G
235//--------------------------
236
237 double c1 = t1[1] ;
238
239 double som = 0;
240 ff0[n3f] = 0 ;
241 for ( i = 3; i < nt; i += 2 ) {
242 ff0[ n3f*i ] = t1[i] - c1 ;
243 som += ff0[ n3f*i ] ;
244 }
245
246// Valeur en psi=0 de la partie antisymetrique de F, F_ :
247 double fmoins0 = nm1s2 * c1 + som ;
248
249// Coef. impairs de G
250// NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
251// donnait exactement les coef. des sinus, ce facteur serait -0.5.
252 for ( i = 3; i < nt; i += 2 ) {
253 g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
254 }
255
256
257// Coefficients pairs de G
258//------------------------
259// Ces coefficients sont egaux aux coefficients pairs du developpement de
260// f.
261// NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
262// donnait exactement les coef. des cosinus, ce facteur serait 1.
263
264 g.set(0) = t1[0] ;
265 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
266 g.set(nm1s2) = t1[nm1] ;
267
268// Transformation de Fourier inverse de G
269//---------------------------------------
270
271// FFT inverse
272 fftw_execute(p) ;
273
274// Valeurs de f deduites de celles de G
275//-------------------------------------
276
277 for ( i = 1; i < nm1s2 ; i++ ) {
278// ... indice du pt symetrique de psi par rapport a pi/2:
279 int isym = nm1 - i ;
280
281 double fp = 0.5 * ( g(i) + g(isym) ) ;
282 double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
283 ff0[ n3f*i ] = ( fp + fm ) / x[isym] ;
284 ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
285 }
286
287//... cas particuliers:
288 ff0[0] = g(0) + fmoins0 ;
289 ff0[ n3f*nm1 ] = 0 ;
290 ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
291
292
293 } // fin de la boucle sur r
294 } // fin de la boucle sur phi
295
296 delete [] t1 ;
297
298}
299}
Lorene prototypes.
Definition app_hor.h:64