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