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 char citsinp_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citsinp.C,v 1.3 2014/10/13 08:53:21 j_novak Exp $" ;
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.3 2014/10/13 08:53:21 j_novak Exp $
88  * $Log: citsinp.C,v $
89  * Revision 1.3 2014/10/13 08:53:21 j_novak
90  * Lorene classes and functions now belong to the namespace Lorene.
91  *
92  * Revision 1.2 2014/10/06 15:18:50 j_novak
93  * Modified #include directives to use c++ syntax.
94  *
95  * Revision 1.1 2004/12/21 17:06:03 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:05 hyc
116  * *** empty log message ***
117  *
118  *
119  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citsinp.C,v 1.3 2014/10/13 08:53:21 j_novak Exp $
120  *
121  */
122 
123 // headers du C
124 #include <cstdlib>
125 #include <fftw3.h>
126 
127 //Lorene prototypes
128 #include "tbl.h"
129 
130 // Prototypage des sous-routines utilisees:
131 namespace Lorene {
132 fftw_plan back_fft(int, Tbl*&) ;
133 double* cheb_ini(const int) ;
134 //*****************************************************************************
135 
136 void citsinp(const int* deg, const int* dimc, double* cf, const int* dimf,
137  double* ff)
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 // Nombres de degres de liberte en theta :
151  int nt = deg[1] ;
152 
153 // Tests de dimension:
154  if (nt > n2f) {
155  cout << "citsinp: nt > n2f : nt = " << nt << " , n2f = "
156  << n2f << endl ;
157  abort () ;
158  exit(-1) ;
159  }
160  if (nt > n2c) {
161  cout << "citsinp: nt > n2c : nt = " << nt << " , n2c = "
162  << n2c << endl ;
163  abort () ;
164  exit(-1) ;
165  }
166  if ( (n1f > 1) && (n1c > n1f) ) {
167  cout << "citsinp: n1c > n1f : n1c = " << n1c << " , n1f = "
168  << n1f << endl ;
169  abort () ;
170  exit(-1) ;
171  }
172  if (n3c > n3f) {
173  cout << "citsinp: n3c > n3f : n3c = " << n3c << " , n3f = "
174  << n3f << 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 = back_fft(nm1, pg) ;
186  Tbl& g = *pg ;
187 
188 // Recherche de la table des sin(psi) :
189  double* sinp = cheb_ini(nt);
190 
191 // boucle sur phi et r (les boucles vont resp. de 0 a max(dimc[0]-2,0)
192 // et 0 a dimc[2]-1)
193 
194  int n2n3f = n2f * n3f ;
195  int n2n3c = n2c * n3c ;
196 
197 /*
198  * Borne de la boucle sur phi:
199  * si n1f = 1, on effectue la boucle une fois seulement.
200  * si n1c > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
201  * j=n1c-1 et j=0 ne sont pas consideres car nuls).
202  */
203  int borne_phi = n1c-1 ;
204  if (n1f == 1) borne_phi = 1 ;
205 
206  for (j=0; j< borne_phi; j++) {
207 
208  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
209 
210  for (k=0; k<n3c; k++) {
211 
212  int i0 = n2n3c * j + k ; // indice de depart
213  double* cf0 = cf + i0 ; // tableau des donnees a transformer
214 
215  i0 = n2n3f * j + k ; // indice de depart
216  double* ff0 = ff + i0 ; // tableau resultat
217 
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 
225 // Calcul des coefficients de Fourier de la fonction
226 // G(psi) = F+(psi) + F_(psi) sin(psi)
227 // en fonction des coefficients en sin(2l theta) de f:
228 
229 //@@
230 // Coefficients
231 //@@
232 // Coefficients en sinus de G
233 //---------------------------
234 // Ces coefficients sont egaux aux coefficients pairs du developmt. en
235 // sin(2l theta) de f (le facteur -.5 vient de la normalisation
236 // de fftw: si fftw donnait reellement les coefficients en sinus,
237 // il faudrait le remplacer par un +1) :
238 
239  for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;
240 
241 // Coefficients en cosinus de G
242 //-----------------------------
243 // Ces coefficients se deduisent des coefficients impairs du developmt.
244 // en sin(2l theta) de f (le facteur +.25 vient de la normalisation
245 // de fftw: si fftw donnait reellement les coefficients en cosinus,
246 // il faudrait le remplacer par un +.5)
247 
248  g.set(0) = .5 * cf0[n3c] ;
249  for ( i = 1; i < nm1s2; i++ ) {
250  g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
251  }
252  g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
253 
254 
255 // Transformation de Fourier inverse de G
256 //---------------------------------------
257 
258 // FFT inverse
259  fftw_execute(p) ;
260 
261 // Valeurs de f deduites de celles de G
262 //-------------------------------------
263 
264  for ( i = 1; i < nm1s2 ; i++ ) {
265 // ... indice du pt symetrique de psi par rapport a pi/2:
266  int isym = nm1 - i ;
267 
268  double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
269  double fm = 0.5 * ( g(i) - g(isym) ) ;
270  ff0[ n3f*i ] = fp + fm ;
271  ff0[ n3f*isym ] = fp - fm ;
272  }
273 
274 //... cas particuliers:
275  ff0[0] = 0. ;
276  ff0[ n3f*nm1 ] = -2. * g(0) ;
277  ff0[ n3f*nm1s2 ] = g(nm1s2) ;
278 
279 
280  } // fin de la boucle sur r
281  } // fin de la boucle sur phi
282 
283 }
284 }
Lorene prototypes.
Definition: app_hor.h:64