LORENE
FFTW3/circhebpip.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 circhebpip_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/circhebpip.C,v 1.3 2014/10/13 08:53:20 j_novak Exp $" ;
24 
25 
26 /*
27  * Transformation de Tchebyshev inverse (cas rare) sur le troisieme indice
28  * (indice correspondant a r) d'un tableau 3-D decrivant une fonction quelconque.
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 r est nr = deg[2] et doit etre de la forme
36  * nr = 2*p + 1
37  * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
38  * dimensions.
39  * On doit avoir dimc[2] >= deg[2] = nr.
40  * NB: pour dimc[0] = 1 (un seul point en phi), la transformation
41  * est bien effectuee.
42  * pour dimc[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 != dimc[0]-1 (cf. commentaires sur borne_phi).
45  *
46  * double* cf : tableau des coefficients c_i de la fonction definis
47  * comme suit (a theta et phi fixes)
48  *
49  * Si l pair f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) ,
50  * Si l impair f(x) = som_{i=0}^{nr-1} c_i T_{2i+1}(x) ,
51  *
52  * ou T_{i}(x) designe le polynome de Tchebyshev de degre i.
53  * Les coefficients c_i (0 <= i <= nr-1) doivent etre stokes
54  * dans le tableau cf comme suit
55  * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ]
56  * ou j et k sont les indices correspondant a phi et theta
57  * respectivement.
58  * L'espace memoire correspondant a ce pointeur doit etre
59  * dimc[0]*dimc[1]*dimc[2] et doit etre alloue avant l'appel a
60  * la routine.
61  *
62  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
63  * dimensions.
64  * On doit avoir dimf[2] >= deg[2] = nr.
65  *
66  * Sortie:
67  * -------
68  * double* ff : tableau des valeurs de la fonction aux nr points de
69  * de collocation
70  *
71  * x_i = sin( pi/2 i/(nr-1) ) 0 <= i <= nr-1
72  *
73  * Les valeurs de la fonction sont stokees dans le
74  * tableau ff comme suit
75  * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ]
76  * ou j et k sont les indices correspondant a phi et theta
77  * respectivement.
78  * L'espace memoire correspondant a ce pointeur doit etre
79  * dimf[0]*dimf[1]*dimf[2] et doit avoir ete alloue avant
80  * l'appel a la routine.
81  *
82  * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
83  * seul tableau, qui constitue une entree/sortie.
84  */
85 
86 /*
87  * $Id: circhebpip.C,v 1.3 2014/10/13 08:53:20 j_novak Exp $
88  * $Log: circhebpip.C,v $
89  * Revision 1.3 2014/10/13 08:53:20 j_novak
90  * Lorene classes and functions now belong to the namespace Lorene.
91  *
92  * Revision 1.2 2014/10/06 15:18:49 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.1 2004/11/23 15:13:50 m_forot
99  * Added the bases for the cases without any equatorial symmetry
100  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
101  *
102  *
103  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/circhebpip.C,v 1.3 2014/10/13 08:53:20 j_novak Exp $
104  *
105  */
106 
107 // headers du C
108 #include <cstdlib>
109 #include <fftw3.h>
110 
111 //Lorene prototypes
112 #include "tbl.h"
113 
114 // Prototypage des sous-routines utilisees:
115 namespace Lorene {
116 fftw_plan back_fft(int, Tbl*&) ;
117 double* cheb_ini(const int) ;
118 double* chebimp_ini(const int ) ;
119 
120 //*****************************************************************************
121 
122 void circhebpip(const int* deg, const int* dimc, double* cf,
123  const int* dimf, double* ff)
124 
125 {
126 
127 int i, j, k ;
128 
129 // Dimensions des tableaux ff et cf :
130  int n1f = dimf[0] ;
131  int n2f = dimf[1] ;
132  int n3f = dimf[2] ;
133  int n1c = dimc[0] ;
134  int n2c = dimc[1] ;
135  int n3c = dimc[2] ;
136 
137 // Nombres de degres de liberte en r :
138  int nr = deg[2] ;
139 
140 // Tests de dimension:
141  if (nr > n3c) {
142  cout << "circhebpip: nr > n3c : nr = " << nr << " , n3c = "
143  << n3c << endl ;
144  abort () ;
145  exit(-1) ;
146  }
147  if (nr > n3f) {
148  cout << "circhebpip: nr > n3f : nr = " << nr << " , n3f = "
149  << n3f << endl ;
150  abort () ;
151  exit(-1) ;
152  }
153  if (n1c > n1f) {
154  cout << "circhebpip: n1c > n1f : n1c = " << n1c << " , n1f = "
155  << n1f << endl ;
156  abort () ;
157  exit(-1) ;
158  }
159  if (n2c > n2f) {
160  cout << "circhebpip: n2c > n2f : n2c = " << n2c << " , n2f = "
161  << n2f << endl ;
162  abort () ;
163  exit(-1) ;
164  }
165 
166 // Nombre de points pour la FFT:
167  int nm1 = nr - 1;
168  int nm1s2 = nm1 / 2;
169 
170 // Recherche des tables pour la FFT:
171  Tbl* pg = 0x0 ;
172  fftw_plan p = back_fft(nm1, pg) ;
173  Tbl& g = *pg ;
174  double* t1 = new double[nr] ;
175 
176 // Recherche de la table des sin(psi) :
177  double* sinp = cheb_ini(nr);
178 
179 // Recherche de la table des points de collocations x_k :
180  double* x = chebimp_ini(nr);
181 
182 // boucle sur phi et theta
183 
184  int n2n3f = n2f * n3f ;
185  int n2n3c = n2c * n3c ;
186 
187 /*
188  * Borne de la boucle sur phi:
189  * si n1c = 1, on effectue la boucle une fois seulement.
190  * si n1c > 1, on va jusqu'a j = n1c-2 en sautant j = 1 (les coefficients
191  * j=n1c-1 et j=0 ne sont pas consideres car nuls).
192  */
193  int borne_phi = ( n1c > 1 ) ? n1c-1 : 1 ;
194 
195  for (j=0; j< borne_phi; j++) {
196 
197  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
198 
199 
200  /************ Cas l pair **********/
201 
202  for (k=0; k<n2c; k+=2) {
203 
204  int i0 = n2n3c * j + n3c * k ; // indice de depart
205  double* cf0 = cf + i0 ; // tableau des donnees a transformer
206 
207  i0 = n2n3f * j + n3f * k ; // indice de depart
208  double* ff0 = ff + i0 ; // tableau resultat
209 
210 
211 // Calcul des coefficients de Fourier de la fonction
212 // G(psi) = F+(theta) + F_(theta) sin(theta)
213 // en fonction des coefficients de Tchebyshev de f:
214 
215 // Coefficients impairs de G
216 //--------------------------
217 
218  double c1 = cf0[1] ;
219 
220  double som = 0;
221  ff0[1] = 0 ;
222  for ( i = 3; i < nr; i += 2 ) {
223  ff0[i] = cf0[i] - c1 ;
224  som += ff0[i] ;
225  }
226 
227 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
228  double fmoins0 = nm1s2 * c1 + som ;
229 
230 // Coef. impairs de G
231 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
232 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
233  for ( i = 3; i < nr; i += 2 ) {
234  g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ;
235  }
236 
237 
238 // Coefficients pairs de G
239 //------------------------
240 // Ces coefficients sont egaux aux coefficients pairs du developpement de
241 // f.
242 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
243 // donnait exactement les coef. des cosinus, ce facteur serait 1.
244 
245  g.set(0) = cf0[0] ;
246  for (i=1; i<nm1s2; i++) g.set(i) = 0.5 * cf0[2*i] ;
247  g.set(nm1s2) = cf0[nm1] ;
248 
249 // Transformation de Fourier inverse de G
250 //---------------------------------------
251 
252 // FFT inverse
253  fftw_execute(p) ;
254 
255 // Valeurs de f deduites de celles de G
256 //-------------------------------------
257 
258  for ( i = 1; i < nm1s2 ; i++ ) {
259 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
260  int isym = nm1 - i ;
261 // ... indice (dans le tableau ff0) du point x correspondant a psi
262  int ix = nm1 - i ;
263 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
264  int ixsym = nm1 - isym ;
265 
266  double fp = .5 * ( g(i) + g(isym) ) ;
267  double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ;
268 
269  ff0[ix] = fp + fm ;
270  ff0[ixsym] = fp - fm ;
271  }
272 
273 //... cas particuliers:
274  ff0[0] = g(0) - fmoins0 ;
275  ff0[nm1] = g(0) + fmoins0 ;
276  ff0[nm1s2] = g(nm1s2) ;
277 
278  } // fin de la boucle sur theta
279 
280 
281  /*********** Cas l impair **********/
282 
283  for (k=1; k<n2c; k+=2) {
284 
285  int i0 = n2n3c * j + n3c * k ; // indice de depart
286  double* cf0 = cf + i0 ; // tableau des donnees a transformer
287 
288  i0 = n2n3f * j + n3f * k ; // indice de depart
289  double* ff0 = ff + i0 ; // tableau resultat
290 
291 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
292 // h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le
293 // tableau t1 :
294  t1[0] = .5 * cf0[0] ;
295  for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ;
296  t1[nm1] = .5 * cf0[nr-2] ;
297 
298 
299 // Calcul des coefficients de Fourier de la fonction
300 // G(psi) = F+(theta) + F_(theta) sin(theta)
301 // en fonction des coefficients de Tchebyshev de f:
302 
303 // Coefficients impairs de G
304 //--------------------------
305 
306  double c1 = t1[1] ;
307 
308  double som = 0;
309  ff0[1] = 0 ;
310  for ( i = 3; i < nr; i += 2 ) {
311  ff0[i] = t1[i] - c1 ;
312  som += ff0[i] ;
313  }
314 
315 // Valeur en theta=0 de la partie antisymetrique de F, F_ :
316  double fmoins0 = nm1s2 * c1 + som ;
317 
318 // Coef. impairs de G
319 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
320 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
321  for ( i = 3; i < nr; i += 2 ) {
322  g.set(nm1-i/2) = 0.25 * ( ff0[i] - ff0[i-2] ) ;
323  }
324 
325 
326 // Coefficients pairs de G
327 //------------------------
328 // Ces coefficients sont egaux aux coefficients pairs du developpement de
329 // f.
330 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
331 // donnait exactement les coef. des cosinus, ce facteur serait 1.
332 
333  g.set(0) = t1[0] ;
334  for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
335  g.set(nm1s2) = t1[nm1] ;
336 
337 // Transformation de Fourier inverse de G
338 //---------------------------------------
339 
340 // FFT inverse
341  fftw_execute(p) ;
342 
343 // Valeurs de f deduites de celles de G
344 //-------------------------------------
345 
346  for ( i = 1; i < nm1s2 ; i++ ) {
347 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
348  int isym = nm1 - i ;
349 // ... indice (dans le tableau ff0) du point x correspondant a psi
350  int ix = nm1 - i ;
351 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
352  int ixsym = nm1 - isym ;
353 
354  double fp = .5 * ( g(i) + g(isym) ) ;
355  double fm = .5 * ( g(i) - g(isym) ) / sinp[i] ;
356 
357  ff0[ix] = ( fp + fm ) / x[ix];
358  ff0[ixsym] = ( fp - fm ) / x[ixsym] ;
359  }
360 
361 //... cas particuliers:
362  ff0[0] = 0 ;
363  ff0[nm1] = g(0) + fmoins0 ;
364  ff0[nm1s2] = g(nm1s2) / x[nm1s2] ;
365 
366  } // fin de la boucle sur theta
367  } // fin de la boucle sur phi
368 
369  delete [] t1 ;
370 
371 }
372 }
373 
Lorene
Lorene prototypes.
Definition: app_hor.h:64