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