LORENE
FFT991/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
23
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 routine FFT Fortran FFT991
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 3^q 5^r + 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.5 2016/12/05 16:18:04 j_novak Exp $
82 * $Log: citsini.C,v $
83 * Revision 1.5 2016/12/05 16:18:04 j_novak
84 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
85 *
86 * Revision 1.4 2014/10/15 12:48:22 j_novak
87 * Corrected namespace declaration.
88 *
89 * Revision 1.3 2014/10/13 08:53:17 j_novak
90 * Lorene classes and functions now belong to the namespace Lorene.
91 *
92 * Revision 1.2 2014/10/06 15:18:47 j_novak
93 * Modified #include directives to use c++ syntax.
94 *
95 * Revision 1.1 2004/12/21 17:06:01 j_novak
96 * Added all files for using fftw3.
97 *
98 * Revision 1.4 2003/01/31 10:31:24 e_gourgoulhon
99 * Suppressed the directive #include <malloc.h> for malloc is defined
100 * in <stdlib.h>
101 *
102 * Revision 1.3 2002/10/16 14:36:54 j_novak
103 * Reorganization of #include instructions of standard C++, in order to
104 * use experimental version 3 of gcc.
105 *
106 * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
107 * Modification of declaration of Fortran 77 prototypes for
108 * a better portability (in particular on IBM AIX systems):
109 * All Fortran subroutine names are now written F77_* and are
110 * defined in the new file C++/Include/proto_f77.h.
111 *
112 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
113 * LORENE
114 *
115 * Revision 2.0 1999/02/22 15:41:16 hyc
116 * *** empty log message ***
117 *
118 *
119 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/citsini.C,v 1.5 2016/12/05 16:18:04 j_novak Exp $
120 *
121 */
122
123// headers du C
124#include <cassert>
125#include <cstdlib>
126
127// Prototypes of F77 subroutines
128#include "headcpp.h"
129#include "proto_f77.h"
130
131// Prototypage des sous-routines utilisees:
132namespace Lorene {
133int* facto_ini(int ) ;
134double* trigo_ini(int ) ;
135double* cheb_ini(const int) ;
136double* chebimp_ini(const int ) ;
137//*****************************************************************************
138
139void citsini(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 << "citsini: nt > n2f : nt = " << nt << " , n2f = "
159 << n2f << endl ;
160 abort () ;
161 exit(-1) ;
162 }
163 if (nt > n2c) {
164 cout << "citsini: nt > n2c : nt = " << nt << " , n2c = "
165 << n2c << endl ;
166 abort () ;
167 exit(-1) ;
168 }
169 if ( (n1f > 1) && (n1c > n1f) ) {
170 cout << "citsini: n1c > n1f : n1c = " << n1c << " , n1f = "
171 << n1f << endl ;
172 abort () ;
173 exit(-1) ;
174 }
175 if (n3c > n3f) {
176 cout << "citsini: 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 int* facto = facto_ini(nm1) ;
188 double* trigo = trigo_ini(nm1) ;
189
190// Recherche de la table des sin(psi) :
191 double* sinp = cheb_ini(nt);
192
193// Recherche de la table des sin( theta_l ) :
194 double* sinth = chebimp_ini(nt);
195
196 // tableau de travail t1 et g
197 // (la dimension nm1+2 = nt+1 est exigee par la routine fft991)
198 double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
199 double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
200
201// Parametres pour la routine FFT991
202 int jump = 1 ;
203 int inc = 1 ;
204 int lot = 1 ;
205 int isign = 1 ;
206
207// boucle sur phi et r
208
209 int n2n3f = n2f * n3f ;
210 int n2n3c = n2c * n3c ;
211
212/*
213 * Borne de la boucle sur phi:
214 * si n1f = 1, on effectue la boucle une fois seulement.
215 * si n1f > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
216 * j=n1c-1 et j=0 ne sont pas consideres car nuls).
217 */
218
219 int borne_phi = n1c-1 ;
220 if (n1f == 1) borne_phi = 1 ;
221
222 for (j=0; j< borne_phi; j++) {
223
224 if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
225
226 for (k=0; k<n3c; k++) {
227
228 int i0 = n2n3c * j + k ; // indice de depart
229 double* cf0 = cf + i0 ; // tableau des donnees a transformer
230
231 i0 = n2n3f * j + k ; // indice de depart
232 double* ff0 = ff + i0 ; // tableau resultat
233
234// Calcul des coefficients du developpement en cos(2 l theta)
235// de la fonction h(theta) := f(theta) sin(theta)
236// en fonction de ceux de f (le resultat est stoke dans le tableau t1) :
237 t1[0] = .5 * cf0[0] ;
238 for (i=1; i<nm1; i++) {
239 t1[i] = .5 * ( cf0[ n3c*i ] - cf0[ n3c*(i-1) ] ) ;
240 }
241 t1[nm1] = -.5 * cf0[ n3c*(nt-2) ] ;
242
243/*
244 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
245 * reliee a theta par psi = 2 theta et F(psi) = h(theta(psi)).
246 */
247
248// Calcul des coefficients de Fourier de la fonction
249// G(psi) = F+(psi) + F_(psi) sin(psi)
250// en fonction des coefficients en cos(2l theta) de h:
251
252// Coefficients impairs de G
253//--------------------------
254
255 double c1 = t1[1] ;
256
257 double som = 0;
258 ff0[n3f] = 0 ;
259 for ( i = 3; i < nt; i += 2 ) {
260 ff0[ n3f*i ] = t1[i] - c1 ;
261 som += ff0[ n3f*i ] ;
262 }
263
264// Valeur en psi=0 de la partie antisymetrique de F, F_ :
265 double fmoins0 = nm1s2 * c1 + som ;
266
267// Coef. impairs de G
268// NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
269// donnait exactement les coef. des sinus, ce facteur serait -0.5.
270 g[1] = 0 ;
271 for ( i = 3; i < nt; i += 2 ) {
272 g[i] = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
273 }
274 g[nt] = 0 ;
275
276
277// Coefficients pairs de G
278//------------------------
279// Ces coefficients sont egaux aux coefficients pairs du developpement de
280// h.
281// NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
282// donnait exactement les coef. des cosinus, ce facteur serait 1.
283
284 g[0] = t1[0] ;
285 for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * t1[i] ;
286 g[nm1] = t1[nm1] ;
287
288// Transformation de Fourier inverse de G
289//---------------------------------------
290
291// FFT inverse
292 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
293
294// Valeurs de f deduites de celles de G
295//-------------------------------------
296
297 for ( i = 1; i < nm1s2 ; i++ ) {
298// ... indice du pt symetrique de psi par rapport a pi/2:
299 int isym = nm1 - i ;
300
301 double fp = 0.5 * ( g[i] + g[isym] ) ;
302 double fm = 0.5 * ( g[i] - g[isym] ) / sinp[i] ;
303 ff0[ n3f*i ] = ( fp + fm ) / sinth[i] ;
304 ff0[ n3f*isym ] = ( fp - fm ) / sinth[isym] ;
305 }
306
307//... cas particuliers:
308 ff0[0] = 0 ;
309 ff0[ n3f*nm1 ] = g[0] - fmoins0 ;
310 ff0[ n3f*nm1s2 ] = g[nm1s2] / sinth[nm1s2];
311
312
313 } // fin de la boucle sur r
314 } // fin de la boucle sur phi
315
316 // Menage
317 free (t1) ;
318 free (g) ;
319
320}
321}
Lorene prototypes.
Definition app_hor.h:67
Coord sinp
Definition map.h:735