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