LORENE
FFTW3/cftcos.C
1 /*
2  * Copyright (c) 1999-2002 Eric Gourgoulhon
3  * Copyright (c) 2002 Jerome Novak
4  *
5  * This file is part of LORENE.
6  *
7  * LORENE is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * LORENE is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with LORENE; if not, write to the Free Software
19  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20  *
21  */
22 
23 /*
24  * Transformation en cos(l*theta) sur le deuxieme indice (theta)
25  * d'un tableau 3-D representant une fonction quelconque (theta
26  * varie entre 0 et pi).Utilise la bibliotheque fftw.
27  *
28  * Entree:
29  * -------
30  * int* deg : tableau du nombre effectif de degres de liberte dans chacune
31  * des 3 dimensions: le nombre de points de collocation
32  * en theta est nt = deg[1] et doit etre de la forme
33  * nt = 2*p + 1
34  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
35  * dimensions.
36  * On doit avoir dimf[1] >= deg[1] = nt.
37  * NB: pour dimf[0] = 1 (un seul point en phi), la transformation
38  * est bien effectuee.
39  * pour dimf[0] > 1 (plus d'un point en phi), la
40  * transformation n'est effectuee que pour les indices (en phi)
41  * j != 1 et j != dimf[0]-1 (cf. commentaires sur borne_phi).
42  *
43  * double* ff : tableau des valeurs de la fonction aux nt points de
44  * de collocation
45  *
46  * theta_l = pi l/(nt-1) 0 <= l <= nt-1
47  *
48  * L'espace memoire correspondant a ce
49  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
50  * etre alloue avant l'appel a la routine.
51  * Les valeurs de la fonction doivent etre stokees
52  * dans le tableau ff comme suit
53  * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
54  * ou j et k sont les indices correspondant a
55  * phi et r respectivement.
56  *
57  * int* dimc : tableau du nombre d'elements de cc dans chacune des trois
58  * dimensions.
59  * On doit avoir dimc[1] >= deg[1] = nt.
60  * Sortie:
61  * -------
62  * double* cf : tableau des coefficients c_l de la fonction definis
63  * comme suit (a r et phi fixes)
64  *
65  * f(theta) = som_{l=0}^{nt-1} c_l cos( l theta ) .
66  *
67  * L'espace memoire correspondant a ce
68  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
69  * etre alloue avant l'appel a la routine.
70  * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
71  * le tableau cf comme suit
72  * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
73  * ou j et k sont les indices correspondant a
74  * phi et r respectivement.
75  *
76  * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
77  * seul tableau, qui constitue une entree/sortie.
78  *
79  */
80 
81 
82 
83 /*
84  * $Id: cftcos.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
85  * $Log: cftcos.C,v $
86  * Revision 1.4 2016/12/05 16:18:05 j_novak
87  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
88  *
89  * Revision 1.3 2014/10/13 08:53:18 j_novak
90  * Lorene classes and functions now belong to the namespace Lorene.
91  *
92  * Revision 1.2 2014/10/06 15:18:48 j_novak
93  * Modified #include directives to use c++ syntax.
94  *
95  * Revision 1.1 2004/12/21 17:06:02 j_novak
96  * Added all files for using fftw3.
97  *
98  * Revision 1.3 2004/10/04 13:42:36 j_novak
99  * Using new and delete instead of malloc.
100  *
101  * Revision 1.2 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.1 2002/11/12 17:43:53 j_novak
106  * Added transformatin functions for T_COS basis.
107  *
108  *
109  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcos.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
110  *
111  */
112 
113 // headers du C
114 #include <cstdlib>
115 #include <fftw3.h>
116 
117 //Lorene prototypes
118 #include "tbl.h"
119 
120 // Prototypage des sous-routines utilisees:
121 namespace Lorene {
122 fftw_plan prepare_fft(int, Tbl*&) ;
123 double* cheb_ini(const int) ;
124 
125 //*****************************************************************************
126 
127 void cftcos(const int* deg, const int* dimf, double* ff, const int* dimc,
128  double* cf)
129 {
130 
131 int i, j, k ;
132 
133 // Dimensions des tableaux ff et cf :
134  int n1f = dimf[0] ;
135  int n2f = dimf[1] ;
136  int n3f = dimf[2] ;
137  int n1c = dimc[0] ;
138  int n2c = dimc[1] ;
139  int n3c = dimc[2] ;
140 
141 // Nombre de degres de liberte en theta :
142  int nt = deg[1] ;
143 
144 // Tests de dimension:
145  if (nt > n2f) {
146  cout << "cftcos: nt > n2f : nt = " << nt << " , n2f = "
147  << n2f << endl ;
148  abort () ;
149  exit(-1) ;
150  }
151  if (nt > n2c) {
152  cout << "cftcos: nt > n2c : nt = " << nt << " , n2c = "
153  << n2c << endl ;
154  abort () ;
155  exit(-1) ;
156  }
157  if (n1f > n1c) {
158  cout << "cftcos: n1f > n1c : n1f = " << n1f << " , n1c = "
159  << n1c << endl ;
160  abort () ;
161  exit(-1) ;
162  }
163  if (n3f > n3c) {
164  cout << "cftcos: n3f > n3c : n3f = " << n3f << " , n3c = "
165  << n3c << endl ;
166  abort () ;
167  exit(-1) ;
168  }
169 
170 // Nombre de points pour la FFT:
171  int nm1 = nt - 1;
172  int nm1s2 = nm1 / 2;
173 
174 // Recherche des tables pour la FFT:
175  Tbl* pg = 0x0 ;
176  fftw_plan p = prepare_fft(nm1, pg) ;
177  Tbl& g = *pg ;
178 
179 // Recherche de la table des sin(psi) :
180  double* sinp = cheb_ini(nt);
181 
182 // boucle sur phi et r (les boucles vont resp. de 0 a max(dimf[0]-2,0) et
183 // 0 a dimf[2]-1 )
184 
185  int n2n3f = n2f * n3f ;
186  int n2n3c = n2c * n3c ;
187 
188 /*
189  * Borne de la boucle sur phi:
190  * si n1f = 1, on effectue la boucle une fois seulement.
191  * si n1f > 1, on va jusqu'a j = n1f-2 en sautant j = 1 (les coefficients
192  * j=n1f-1 et j=0 ne sont pas consideres car nuls).
193  */
194  int borne_phi = ( n1f > 1 ) ? n1f-1 : 1 ;
195 
196  for (j=0; j< borne_phi; j++) {
197 
198  if (j==1) continue ; // on ne traite pas le terme en sin(0 phi)
199 
200  for (k=0; k<n3f; k++) {
201 
202  int i0 = n2n3f * j + k ; // indice de depart
203  double* ff0 = ff + i0 ; // tableau des donnees a transformer
204 
205  i0 = n2n3c * j + k ; // indice de depart
206  double* cf0 = cf + i0 ; // tableau resultat
207 
208 /*
209  * NB: dans les commentaires qui suivent, psi designe theta dans [0, pi]
210  */
211 
212 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
213  double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
214 
215 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
216 //---------------------------------------------
217  for ( i = 1; i < nm1s2 ; i++ ) {
218 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
219  int isym = nm1 - i ;
220 // ... indice (dans le tableau ff0) du point theta correspondant a psi
221  int ix = n3f * i ;
222 // ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
223  int ixsym = n3f * isym ;
224 // ... F+(psi)
225  double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
226 // ... F_(psi) sin(psi)
227  double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
228  g.set(i) = fp + fms ;
229  g.set(isym) = fp - fms ;
230  }
231 //... cas particuliers:
232  g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
233  g.set(nm1s2) = ff0[ n3f*nm1s2 ];
234 
235 // Developpement de G en series de Fourier par une FFT
236 //----------------------------------------------------
237 
238  fftw_execute(p) ;
239 
240 // Coefficients pairs du developmt. cos(l theta) de f
241 //----------------------------------------------------
242 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
243 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
244 // de fftw) :
245 
246  double fac = 2./double(nm1) ;
247  cf0[0] = g(0) / double(nm1) ;
248  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
249  cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
250 
251 // Coefficients impairs du developmt. en cos(l theta) de f
252 //---------------------------------------------------------
253 // 1. Coef. c'_k (recurrence amorcee a partir de zero):
254 // Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
255 // (si fftw donnait reellement les coef. en sinus, il faudrait le
256 // remplacer par un -2.)
257  fac *= 2. ;
258  cf0[n3c] = 0 ;
259  double som = 0;
260  for ( i = 3; i < nt; i += 2 ) {
261  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
262  som += cf0[n3c*i] ;
263  }
264 
265 // 2. Calcul de c_1 :
266  double c1 = ( fmoins0 - som ) / nm1s2 ;
267 
268 // 3. Coef. c_k avec k impair:
269  cf0[n3c] = c1 ;
270  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
271 
272 
273  } // fin de la boucle sur r
274  } // fin de la boucle sur phi
275 
276 }
277 }
278 
Lorene prototypes.
Definition: app_hor.h:67