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