LORENE
FFTW3/cftcossins.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
23char cftcossins_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossins.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $" ;
24
25/*
26 * Transformation en sin(l*theta) ou cos(l*theta) (suivant la parite
27 * de l'indice m en phi) sur le deuxieme indice (theta)
28 * d'un tableau 3-D representant une fonction sans symetrie 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 *
42 * double* ff : tableau des valeurs de la fonction aux nt points de
43 * de collocation
44 *
45 * theta_l = pi l/(nt-1) 0 <= l <= nt-1
46 *
47 * L'espace memoire correspondant a ce
48 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
49 * etre alloue avant l'appel a la routine.
50 * Les valeurs de la fonction doivent etre stokees
51 * dans le tableau ff comme suit
52 * f( theta_l ) = ff[ dimf[1]*dimf[2] * m + k + dimf[2] * l ]
53 * ou m et k sont les indices correspondant a
54 * phi et r respectivement.
55 * NB: cette routine suppose que la transformation en phi a deja ete
56 * effectuee: ainsi m est un indice de Fourier, non un indice de
57 * point de collocation en phi.
58 *
59 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
60 * dimensions.
61 * On doit avoir dimc[1] >= deg[1] = nt.
62 * Sortie:
63 * -------
64 * double* cf : tableau des coefficients c_l de la fonction definis
65 * comme suit (a r et phi fixes)
66 *
67 * pour m pair:
68 * f(theta) = som_{l=0}^{nt-1} c_l sin(l theta ) .
69 * pour m impair:
70 * f(theta) = som_{l=0}^{nt-1} c_l cos( l theta ) .
71 *
72 * L'espace memoire correspondant a ce
73 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
74 * etre alloue avant l'appel a la routine.
75 * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
76 * le tableau cf comme suit
77 * c_l = cf[ dimc[1]*dimc[2] * m + k + dimc[2] * l ]
78 * ou m et k sont les indices correspondant a
79 * phi et r respectivement.
80 * Pour m pair, c_0 = c_{nt-1} = 0.
81 * Pour m impair, c_{nt-1} = 0.
82 *
83 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
84 * seul tableau, qui constitue une entree/sortie.
85 *
86 */
87
88/*
89 * $Id: cftcossins.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $
90 * $Log: cftcossins.C,v $
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:48 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.1 2004/11/23 15:13:50 m_forot
101 * Added the bases for the cases without any equatorial symmetry
102 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
103 *
104 *
105 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossins.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $
106 *
107 */
108
109// headers du C
110#include <cstdlib>
111#include <fftw3.h>
112
114#include "tbl.h"
115
116// Prototypage des sous-routines utilisees:
117namespace Lorene {
118fftw_plan prepare_fft(int, Tbl*&) ;
119double* cheb_ini(const int) ;
120//*****************************************************************************
121
122void cftcossins(const int* deg, const int* dimf, double* ff, const int* dimc,
123 double* cf)
124{
125
126int i, j, k ;
127
128// Dimensions des tableaux ff et cf :
129 int n1f = dimf[0] ;
130 int n2f = dimf[1] ;
131 int n3f = dimf[2] ;
132 int n1c = dimc[0] ;
133 int n2c = dimc[1] ;
134 int n3c = dimc[2] ;
135
136// Nombre de degres de liberte en theta :
137 int nt = deg[1] ;
138
139// Tests de dimension:
140 if (nt > n2f) {
141 cout << "cftcossins: nt > n2f : nt = " << nt << " , n2f = "
142 << n2f << endl ;
143 abort () ;
144 exit(-1) ;
145 }
146 if (nt > n2c) {
147 cout << "cftcossins: nt > n2c : nt = " << nt << " , n2c = "
148 << n2c << endl ;
149 abort () ;
150 exit(-1) ;
151 }
152 if (n1f > n1c) {
153 cout << "cftcossins: n1f > n1c : n1f = " << n1f << " , n1c = "
154 << n1c << endl ;
155 abort () ;
156 exit(-1) ;
157 }
158 if (n3f > n3c) {
159 cout << "cftcossins: n3f > n3c : n3f = " << n3f << " , n3c = "
160 << n3c << endl ;
161 abort () ;
162 exit(-1) ;
163 }
164
165// Nombre de points pour la FFT:
166 int nm1 = nt - 1;
167 int nm1s2 = nm1 / 2;
168
169// Recherche des tables pour la FFT:
170 Tbl* pg = 0x0 ;
171 fftw_plan p = prepare_fft(nm1, pg) ;
172 Tbl& g = *pg ;
173
174// Recherche de la table des sin(psi) :
175 double* sinp = cheb_ini(nt);
176
177// boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
178// et 0 a dimf[2])
179
180 int n2n3f = n2f * n3f ;
181 int n2n3c = n2c * n3c ;
182
183//=======================================================================
184// Cas m pair
185//=======================================================================
186
187 j = 0 ;
188
189 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
190 // (car nul)
191
192//--------------------------------------------------------------------
193// partie cos(m phi) avec m pair : transformation en sin(l) theta)
194//--------------------------------------------------------------------
195
196 for (k=0; k<n3f; k++) {
197
198 int i0 = n2n3f * j + k ; // indice de depart
199 double* ff0 = ff + i0 ; // tableau des donnees a transformer
200
201 i0 = n2n3c * j + k ; // indice de depart
202 double* cf0 = cf + i0 ; // tableau resultat
203
204// Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
205//---------------------------------------------
206
207 for ( i = 1; i < nm1s2 ; i++ ) {
208 int isym = nm1 - i ;
209 int ix = n3f * i ;
210 int ixsym = n3f * isym ;
211// ... F+(theta)
212 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
213// ... F_(theta) sin(theta)
214 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
215 g.set(i) = fp + fms ;
216 g.set(isym) = fp - fms ;
217 }
218//... cas particuliers:
219 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
220 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
221
222// Developpement de G en series de Fourier par une FFT
223//----------------------------------------------------
224
225 fftw_execute(p) ;
226
227// Coefficients pairs du developmt. sin(l theta) de f
228//----------------------------------------------------
229// Ces coefficients sont egaux aux coefficients en sinus du developpement
230// de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
231// de fftw) :
232
233 double fac = -2. / double(nm1) ;
234 cf0[0] = 0. ;
235 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
236 cf0[n3c*nm1] = 0. ;
237
238// Coefficients impairs du developmt. en sin(l theta) de f
239//---------------------------------------------------------
240// 1. Coef. c'_k (recurrence amorcee a partir de zero):
241// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
242// (si fftw donnait reellement les coef. en sinus, il faudrait le
243// remplacer par un -2.)
244
245 cf0[n3c] = -fac * g(0);
246 fac *= -2. ;
247 for ( i = 3; i < nt; i += 2 ) {
248 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
249 }
250
251 } // fin de la boucle sur r
252
253//------------------------------------------------------------------------
254// partie sin(m phi) avec m pair : transformation en sin(l theta)
255//------------------------------------------------------------------------
256
257 j++ ;
258
259 if ( j != n1f-1 ) {
260// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
261// pas nuls
262
263 for (k=0; k<n3f; k++) {
264
265 int i0 = n2n3f * j + k ; // indice de depart
266 double* ff0 = ff + i0 ; // tableau des donnees a transformer
267
268 i0 = n2n3c * j + k ; // indice de depart
269 double* cf0 = cf + i0 ; // tableau resultat
270
271// Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
272//---------------------------------------------
273 for ( i = 1; i < nm1s2 ; i++ ) {
274 int isym = nm1 - i ;
275 int ix = n3f * i ;
276 int ixsym = n3f * isym ;
277// ... F+(theta)
278 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
279// ... F_(theta) sin(theta)
280 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
281 g.set(i) = fp + fms ;
282 g.set(isym) = fp - fms ;
283 }
284//... cas particuliers:
285 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
286 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
287
288// Developpement de G en series de Fourier par une FFT
289//----------------------------------------------------
290
291 fftw_execute(p) ;
292
293// Coefficients pairs du developmt. sin(l theta) de f
294//----------------------------------------------------
295// Ces coefficients sont egaux aux coefficients en sinus du developpement
296// de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
297// de fftw) :
298
299 double fac = -2. / double(nm1) ;
300 cf0[0] = 0. ;
301 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
302 cf0[n3c*nm1] = 0. ;
303
304// Coefficients impairs du developmt. en sin(l theta) de f
305//---------------------------------------------------------
306// 1. Coef. c'_k (recurrence amorcee a partir de zero):
307// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
308// (si fftw donnait reellement les coef. en sinus, il faudrait le
309// remplacer par un -2.)
310
311 cf0[n3c] = -fac * g(0);
312 fac *= -2. ;
313 for ( i = 3; i < nt; i += 2 ) {
314 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
315 }
316
317 } // fin de la boucle sur r
318
319 } // fin du cas ou le calcul etait necessaire (i.e. ou les
320 // coef en phi n'etaient pas nuls)
321
322
323// On passe au cas m impair suivant:
324 j+=3 ;
325
326 } // fin de la boucle sur les cas m pair
327
328 if (n1f<=3) // cas m=0 seulement (symetrie axiale)
329 return ;
330
331//=======================================================================
332// Cas m impair
333//=======================================================================
334
335 j = 2 ;
336
337 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
338 // (car nul)
339
340//------------------------------------------------------------------------
341// partie cos(m phi) avec m impair : transformation en cos(l theta)
342//------------------------------------------------------------------------
343
344
345 for (k=0; k<n3f; k++) {
346
347 int i0 = n2n3f * j + k ; // indice de depart
348 double* ff0 = ff + i0 ; // tableau des donnees a transformer
349
350 i0 = n2n3c * j + k ; // indice de depart
351 double* cf0 = cf + i0 ; // tableau resultat
352
353
354// Valeur en theta=0 de la partie antisymetrique de F, F_ :
355 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
356
357// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
358//---------------------------------------------
359 for ( i = 1; i < nm1s2 ; i++ ) {
360 int isym = nm1 - i ;
361 int ix = n3f * i ;
362 int ixsym = n3f * isym ;
363// ... F+(theta)
364 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
365// ... F_(theta) sin(psi)
366 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
367 g.set(i) = fp + fms ;
368 g.set(isym) = fp - fms ;
369 }
370//... cas particuliers:
371 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
372 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
373
374// Developpement de G en series de Fourier par une FFT
375//----------------------------------------------------
376
377 fftw_execute(p) ;
378
379// Coefficients pairs du developmt. cos(l theta) de f
380//----------------------------------------------------
381// Ces coefficients sont egaux aux coefficients en cosinus du developpement
382// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
383// de fftw) :
384
385 double fac = 2./double(nm1) ;
386 cf0[0] = g(0) / double(nm1) ;
387 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
388 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
389
390// Coefficients impairs du developmt. en cos(l theta) de f
391//---------------------------------------------------------
392// 1. Coef. c'_k (recurrence amorcee a partir de zero):
393// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
394// (si fftw donnait reellement les coef. en sinus, il faudrait le
395// remplacer par un -2.)
396 fac *= 2. ;
397 cf0[n3c] = 0 ;
398 double som = 0;
399 for ( i = 3; i < nt; i += 2 ) {
400 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
401 som += cf0[n3c*i] ;
402 }
403
404// 2. Calcul de c_1 :
405 double c1 = ( fmoins0 - som ) / nm1s2 ;
406
407// 3. Coef. c_k avec k impair:
408 cf0[n3c] = c1 ;
409 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
410
411
412 } // fin de la boucle sur r
413
414//--------------------------------------------------------------------
415// partie sin(m phi) avec m impair : transformation en cos(l theta)
416//--------------------------------------------------------------------
417
418 j++ ;
419
420 if ( (j != 1) && (j != n1f-1 ) ) {
421// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
422// pas nuls
423
424 for (k=0; k<n3f; k++) {
425
426 int i0 = n2n3f * j + k ; // indice de depart
427 double* ff0 = ff + i0 ; // tableau des donnees a transformer
428
429 i0 = n2n3c * j + k ; // indice de depart
430 double* cf0 = cf + i0 ; // tableau resultat
431
432
433// Valeur en theta=0 de la partie antisymetrique de F, F_ :
434 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
435
436// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
437//---------------------------------------------
438 for ( i = 1; i < nm1s2 ; i++ ) {
439 int isym = nm1 - i ;
440 int ix = n3f * i ;
441 int ixsym = n3f * isym ;
442// ... F+(theta)
443 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
444// ... F_(theta) sin(psi)
445 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
446 g.set(i) = fp + fms ;
447 g.set(isym) = fp - fms ;
448 }
449//... cas particuliers:
450 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
451 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
452
453// Developpement de G en series de Fourier par une FFT
454//----------------------------------------------------
455
456 fftw_execute(p) ;
457
458// Coefficients pairs du developmt. cos(l theta) de f
459//----------------------------------------------------
460// Ces coefficients sont egaux aux coefficients en cosinus du developpement
461// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
462// de fftw) :
463
464 double fac = 2./double(nm1) ;
465 cf0[0] = g(0) / double(nm1) ;
466 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
467 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
468
469// Coefficients impairs du developmt. en cos(l theta) de f
470//---------------------------------------------------------
471// 1. Coef. c'_k (recurrence amorcee a partir de zero):
472// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
473// (si fftw donnait reellement les coef. en sinus, il faudrait le
474// remplacer par un -2.)
475 fac *= 2. ;
476 cf0[n3c] = 0 ;
477 double som = 0;
478 for ( i = 3; i < nt; i += 2 ) {
479 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
480 som += cf0[n3c*i] ;
481 }
482
483// 2. Calcul de c_1 :
484 double c1 = ( fmoins0 - som ) / nm1s2 ;
485
486// 3. Coef. c_k avec k impair:
487 cf0[n3c] = c1 ;
488 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
489
490
491 } // fin de la boucle sur r
492 } // fin du cas ou le calcul etait necessaire (i.e. ou les
493 // coef en phi n'etaient pas nuls)
494
495// On passe au cas m impair suivant:
496 j+=3 ;
497
498 } // fin de la boucle sur les cas m pair
499
500}
501}
Lorene prototypes.
Definition app_hor.h:64