LORENE
FFT991/cftcosi.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 cftcosi_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cftcosi.C,v 1.4 2014/10/15 12:48:20 j_novak Exp $" ;
24
25/*
26 * Transformation en cos((2*l+1)*theta) 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 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* dimf : tableau du nombre d'elements de ff dans chacune des trois
38 * dimensions.
39 * On doit avoir dimf[1] >= deg[1] = nt.
40 * NB: pour dimf[0] = 1 (un seul point en phi), la transformation
41 * est bien effectuee.
42 * pour dimf[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 != dimf[0]-1 (cf. commentaires sur borne_phi).
45 *
46 * double* ff : tableau des valeurs de la fonction aux nt points de
47 * de collocation
48 *
49 * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
50 *
51 * L'espace memoire correspondant a ce
52 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
53 * etre alloue avant l'appel a la routine.
54 * Les valeurs de la fonction doivent etre stokees
55 * dans le tableau ff comme suit
56 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
57 * ou j et k sont les indices correspondant a
58 * phi et r respectivement.
59 *
60 * int* dimc : tableau du nombre d'elements de cc dans chacune des trois
61 * dimensions.
62 * On doit avoir dimc[1] >= deg[1] = nt.
63 * Sortie:
64 * -------
65 * double* cf : tableau des coefficients c_l de la fonction definis
66 * comme suit (a r et phi fixes)
67 *
68 * f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) .
69 *
70 * L'espace memoire correspondant a ce
71 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
72 * etre alloue avant l'appel a la routine.
73 * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
74 * le tableau cf comme suit
75 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
76 * ou j et k sont les indices correspondant a
77 * phi et r respectivement. On a c_{nt-1}=0.
78 *
79 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
80 * seul tableau, qui constitue une entree/sortie.
81 *
82 */
83
84/*
85 * $Id: cftcosi.C,v 1.4 2014/10/15 12:48:20 j_novak Exp $
86 * $Log: cftcosi.C,v $
87 * Revision 1.4 2014/10/15 12:48:20 j_novak
88 * Corrected namespace declaration.
89 *
90 * Revision 1.3 2014/10/13 08:53:15 j_novak
91 * Lorene classes and functions now belong to the namespace Lorene.
92 *
93 * Revision 1.2 2014/10/06 15:18:45 j_novak
94 * Modified #include directives to use c++ syntax.
95 *
96 * Revision 1.1 2004/12/21 17:06:01 j_novak
97 * Added all files for using fftw3.
98 *
99 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
100 * Suppressed the directive #include <malloc.h> for malloc is defined
101 * in <stdlib.h>
102 *
103 * Revision 1.3 2002/10/16 14:36:44 j_novak
104 * Reorganization of #include instructions of standard C++, in order to
105 * use experimental version 3 of gcc.
106 *
107 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
108 * Modification of declaration of Fortran 77 prototypes for
109 * a better portability (in particular on IBM AIX systems):
110 * All Fortran subroutine names are now written F77_* and are
111 * defined in the new file C++/Include/proto_f77.h.
112 *
113 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
114 * LORENE
115 *
116 * Revision 2.0 1999/02/22 15:47:57 hyc
117 * *** empty log message ***
118 *
119 *
120 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cftcosi.C,v 1.4 2014/10/15 12:48:20 j_novak Exp $
121 *
122 */
123
124
125// headers du C
126#include <cassert>
127#include <cstdlib>
128
129// Prototypes of F77 subroutines
130#include "headcpp.h"
131#include "proto_f77.h"
132
133// Prototypage des sous-routines utilisees:
134namespace Lorene {
135int* facto_ini(int ) ;
136double* trigo_ini(int ) ;
137double* cheb_ini(const int) ;
138double* chebimp_ini(const int ) ;
139//*****************************************************************************
140
141void cftcosi(const int* deg, const int* dimf, double* ff, const int* dimc,
142 double* cf)
143{
144
145int i, j, k ;
146
147// Dimensions des tableaux ff et cf :
148 int n1f = dimf[0] ;
149 int n2f = dimf[1] ;
150 int n3f = dimf[2] ;
151 int n1c = dimc[0] ;
152 int n2c = dimc[1] ;
153 int n3c = dimc[2] ;
154
155// Nombre de degres de liberte en theta :
156 int nt = deg[1] ;
157
158// Tests de dimension:
159 if (nt > n2f) {
160 cout << "cftcosi: nt > n2f : nt = " << nt << " , n2f = "
161 << n2f << endl ;
162 abort () ;
163 exit(-1) ;
164 }
165 if (nt > n2c) {
166 cout << "cftcosi: nt > n2c : nt = " << nt << " , n2c = "
167 << n2c << endl ;
168 abort () ;
169 exit(-1) ;
170 }
171 if (n1f > n1c) {
172 cout << "cftcosi: n1f > n1c : n1f = " << n1f << " , n1c = "
173 << n1c << endl ;
174 abort () ;
175 exit(-1) ;
176 }
177 if (n3f > n3c) {
178 cout << "cftcosi: n3f > n3c : n3f = " << n3f << " , n3c = "
179 << n3c << endl ;
180 abort () ;
181 exit(-1) ;
182 }
183
184// Nombre de points pour la FFT:
185 int nm1 = nt - 1;
186 int nm1s2 = nm1 / 2;
187
188// Recherche des tables pour la FFT:
189 int* facto = facto_ini(nm1) ;
190 double* trigo = trigo_ini(nm1) ;
191
192// Recherche de la table des sin(psi) :
193 double* sinp = cheb_ini(nt);
194
195// Recherche de la table des points de collocations x_k :
196 double* x = chebimp_ini(nt);
197
198 // tableau de travail G et t1
199 double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) );
200 double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
201
202// Parametres pour la routine FFT991
203 int jump = 1 ;
204 int inc = 1 ;
205 int lot = 1 ;
206 int isign = - 1 ;
207
208// boucle sur phi et r (les boucles vont resp. de 0 a max(dimf[0]-2,0) et
209// 0 a dimf[2]-1 )
210
211 int n2n3f = n2f * n3f ;
212 int n2n3c = n2c * n3c ;
213
214/*
215 * Borne de la boucle sur phi:
216 * si n1f = 1, on effectue la boucle une fois seulement.
217 * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
218 * j=n1f-1 et j=0 ne sont pas consideres car nuls).
219 */
220 int borne_phi = ( n1f > 1 ) ? n1f-1 : 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<n3f; k++) {
227
228 int i0 = n2n3f * j + k ; // indice de depart
229 double* ff0 = ff + i0 ; // tableau des donnees a transformer
230
231 i0 = n2n3c * j + k ; // indice de depart
232 double* cf0 = cf + i0 ; // tableau resultat
233
234// Multiplication de la fonction par x=cos(theta) (pour la rendre paire)
235// NB: dans les commentaires qui suivent, on note h(x) = x f(x).
236// (Les valeurs de h dans l'espace des configurations sont stokees dans le
237// tableau cf0).
238 for (i=0; i<nt-1; i++) cf0[n3c*i] = x[nm1-i] * ff0[n3f*i] ;
239 cf0[n3c*nm1] = 0 ;
240
241/*
242 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
243 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
244 */
245
246// Valeur en psi=0 de la partie antisymetrique de F, F_ :
247 double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
248
249// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
250//---------------------------------------------
251 for ( i = 1; i < nm1s2 ; i++ ) {
252// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
253 int isym = nm1 - i ;
254// ... indice (dans le tableau cf0) du point theta correspondant a psi
255 int ix = n3c * i ;
256// ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
257 int ixsym = n3c * isym ;
258// ... F+(psi)
259 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
260// ... F_(psi) sin(psi)
261 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
262 g[i] = fp + fms ;
263 g[isym] = fp - fms ;
264 }
265//... cas particuliers:
266 g[0] = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
267 g[nm1s2] = cf0[ n3c*nm1s2 ];
268
269// Developpement de G en series de Fourier par une FFT
270//----------------------------------------------------
271
272 F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
273
274// Coefficients pairs du developmt. de Tchebyshev de h
275//----------------------------------------------------
276// Ces coefficients sont egaux aux coefficients en cosinus du developpement
277// de G en series de Fourier (le facteur 2. vient de la normalisation
278// de fft991; si fft991 donnait reellement les coef. en cosinus, il faudrait le
279// remplacer par un +1.) :
280
281 cf0[0] = g[0] ;
282 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = 2.* g[i] ;
283 cf0[n3c*nm1] = g[nm1] ;
284
285// Coefficients impairs du developmt. de Tchebyshev de h
286//------------------------------------------------------
287// 1. Coef. c'_k (recurrence amorcee a partir de zero)
288// Le +4. en facteur de g[i] est du a la normalisation de fft991
289// (si fft991 donnait reellement les coef. en sinus, il faudrait le
290// remplacer par un -2.)
291 cf0[n3c] = 0 ;
292 double som = 0;
293 for ( i = 3; i < nt; i += 2 ) {
294 cf0[n3c*i] = cf0[n3c*(i-2)] + 4. * g[i] ;
295 som += cf0[n3c*i] ;
296 }
297
298// 2. Calcul de c_1 :
299 double c1 = ( fmoins0 - som ) / nm1s2 ;
300
301// 3. Coef. c_k avec k impair:
302 cf0[n3c] = c1 ;
303 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
304
305
306// Coefficients de f en fonction de ceux de h
307//-------------------------------------------
308
309 cf0[0] = 2* cf0[0] ;
310 for (i=1; i<nm1; i++) {
311 cf0[n3c*i] = 2 * cf0[n3c*i] - cf0[n3c*(i-1)] ;
312 }
313 cf0[n3c*nm1] = 0 ;
314
315 } // fin de la boucle sur r
316 } // fin de la boucle sur phi
317
318 // Menage
319 free (t1) ;
320 free (g) ;
321
322}
323}
Lorene prototypes.
Definition app_hor.h:64