LORENE
FFTW3/citsinp.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 sin(2*l*theta) inverse sur le deuxieme indice (theta)
28 * d'un tableau 3-D representant une fonction antisymetrique 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* dimc : tableau du nombre d'elements de cc dans chacune des trois
39 * dimensions.
40 * On doit avoir dimc[1] >= deg[1] = nt.
41 * NB: pour dimc[0] = 1 (un seul point en phi), la transformation
42 * est bien effectuee.
43 * pour dimc[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 != dimc[0]-1 (cf. commentaires sur borne_phi).
46 *
47 * double* cf : tableau des coefficients c_l de la fonction definis
48 * comme suit (a r et phi fixes)
49 *
50 * f(theta) = som_{l=0}^{nt-1} c_l sin( 2 l theta ) .
51 *
52 * L'espace memoire correspondant a ce
53 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
54 * avoir ete alloue avant l'appel a la routine.
55 * Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
56 * le tableau cf comme suit
57 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
58 * ou j et k sont les indices correspondant a
59 * phi et r respectivement.
60 *
61 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
62 * dimensions.
63 * On doit avoir dimf[1] >= deg[1] = nt.
64 *
65 * Sortie:
66 * -------
67 * double* ff : tableau des valeurs de la fonction aux nt points de
68 * de collocation
69 *
70 * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
71 *
72 * L'espace memoire correspondant a ce
73 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
74 * avoir ete alloue avant l'appel a la routine.
75 * Les valeurs de la fonction sont stokees
76 * dans le tableau ff comme suit
77 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
78 * ou j et k sont les indices correspondant a
79 * phi et r respectivement.
80 *
81 * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
82 * seul tableau, qui constitue une entree/sortie.
83 *
84 */
85
86/*
87 * $Id: citsinp.C,v 1.4 2016/12/05 16:18:06 j_novak Exp $
88 * $Log: citsinp.C,v $
89 * Revision 1.4 2016/12/05 16:18:06 j_novak
90 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
91 *
92 * Revision 1.3 2014/10/13 08:53:21 j_novak
93 * Lorene classes and functions now belong to the namespace Lorene.
94 *
95 * Revision 1.2 2014/10/06 15:18:50 j_novak
96 * Modified #include directives to use c++ syntax.
97 *
98 * Revision 1.1 2004/12/21 17:06:03 j_novak
99 * Added all files for using fftw3.
100 *
101 * Revision 1.4 2003/01/31 10:31:24 e_gourgoulhon
102 * Suppressed the directive #include <malloc.h> for malloc is defined
103 * in <stdlib.h>
104 *
105 * Revision 1.3 2002/10/16 14:36:54 j_novak
106 * Reorganization of #include instructions of standard C++, in order to
107 * use experimental version 3 of gcc.
108 *
109 * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
110 * Modification of declaration of Fortran 77 prototypes for
111 * a better portability (in particular on IBM AIX systems):
112 * All Fortran subroutine names are now written F77_* and are
113 * defined in the new file C++/Include/proto_f77.h.
114 *
115 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
116 * LORENE
117 *
118 * Revision 2.0 1999/02/22 15:41:05 hyc
119 * *** empty log message ***
120 *
121 *
122 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citsinp.C,v 1.4 2016/12/05 16:18:06 j_novak Exp $
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) ;
137//*****************************************************************************
138
139void citsinp(const int* deg, const int* dimc, double* cf, const int* dimf,
140 double* ff)
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 theta :
154 int nt = deg[1] ;
155
156// Tests de dimension:
157 if (nt > n2f) {
158 cout << "citsinp: nt > n2f : nt = " << nt << " , n2f = "
159 << n2f << endl ;
160 abort () ;
161 exit(-1) ;
162 }
163 if (nt > n2c) {
164 cout << "citsinp: nt > n2c : nt = " << nt << " , n2c = "
165 << n2c << endl ;
166 abort () ;
167 exit(-1) ;
168 }
169 if ( (n1f > 1) && (n1c > n1f) ) {
170 cout << "citsinp: n1c > n1f : n1c = " << n1c << " , n1f = "
171 << n1f << endl ;
172 abort () ;
173 exit(-1) ;
174 }
175 if (n3c > n3f) {
176 cout << "citsinp: n3c > n3f : n3c = " << n3c << " , n3f = "
177 << n3f << 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 = back_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(dimc[0]-2,0)
195// et 0 a dimc[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 n1c > 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
222/*
223 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
224 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
225 */
226
227
228// Calcul des coefficients de Fourier de la fonction
229// G(psi) = F+(psi) + F_(psi) sin(psi)
230// en fonction des coefficients en sin(2l theta) de f:
231
232//@@
233// Coefficients
234//@@
235// Coefficients en sinus de G
236//---------------------------
237// Ces coefficients sont egaux aux coefficients pairs du developmt. en
238// sin(2l theta) de f (le facteur -.5 vient de la normalisation
239// de fftw: si fftw donnait reellement les coefficients en sinus,
240// il faudrait le remplacer par un +1) :
241
242 for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;
243
244// Coefficients en cosinus de G
245//-----------------------------
246// Ces coefficients se deduisent des coefficients impairs du developmt.
247// en sin(2l theta) de f (le facteur +.25 vient de la normalisation
248// de fftw: si fftw donnait reellement les coefficients en cosinus,
249// il faudrait le remplacer par un +.5)
250
251 g.set(0) = .5 * cf0[n3c] ;
252 for ( i = 1; i < nm1s2; i++ ) {
253 g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
254 }
255 g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
256
257
258// Transformation de Fourier inverse de G
259//---------------------------------------
260
261// FFT inverse
262 fftw_execute(p) ;
263
264// Valeurs de f deduites de celles de G
265//-------------------------------------
266
267 for ( i = 1; i < nm1s2 ; i++ ) {
268// ... indice du pt symetrique de psi par rapport a pi/2:
269 int isym = nm1 - i ;
270
271 double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
272 double fm = 0.5 * ( g(i) - g(isym) ) ;
273 ff0[ n3f*i ] = fp + fm ;
274 ff0[ n3f*isym ] = fp - fm ;
275 }
276
277//... cas particuliers:
278 ff0[0] = 0. ;
279 ff0[ n3f*nm1 ] = -2. * g(0) ;
280 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
281
282
283 } // fin de la boucle sur r
284 } // fin de la boucle sur phi
285
286}
287}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord sinp
Definition map.h:735