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
23
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.4 2016/12/05 16:18:05 j_novak Exp $
90 * $Log: cftcossins.C,v $
91 * Revision 1.4 2016/12/05 16:18:05 j_novak
92 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
93 *
94 * Revision 1.3 2014/10/13 08:53:19 j_novak
95 * Lorene classes and functions now belong to the namespace Lorene.
96 *
97 * Revision 1.2 2014/10/06 15:18:48 j_novak
98 * Modified #include directives to use c++ syntax.
99 *
100 * Revision 1.1 2004/12/21 17:06:02 j_novak
101 * Added all files for using fftw3.
102 *
103 * Revision 1.1 2004/11/23 15:13:50 m_forot
104 * Added the bases for the cases without any equatorial symmetry
105 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
106 *
107 *
108 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossins.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
109 *
110 */
111
112// headers du C
113#include <cstdlib>
114#include <fftw3.h>
115
117#include "tbl.h"
118
119// Prototypage des sous-routines utilisees:
120namespace Lorene {
121fftw_plan prepare_fft(int, Tbl*&) ;
122double* cheb_ini(const int) ;
123//*****************************************************************************
124
125void cftcossins(const int* deg, const int* dimf, double* ff, const int* dimc,
126 double* cf)
127{
128
129int i, j, k ;
130
131// Dimensions des tableaux ff et cf :
132 int n1f = dimf[0] ;
133 int n2f = dimf[1] ;
134 int n3f = dimf[2] ;
135 int n1c = dimc[0] ;
136 int n2c = dimc[1] ;
137 int n3c = dimc[2] ;
138
139// Nombre de degres de liberte en theta :
140 int nt = deg[1] ;
141
142// Tests de dimension:
143 if (nt > n2f) {
144 cout << "cftcossins: nt > n2f : nt = " << nt << " , n2f = "
145 << n2f << endl ;
146 abort () ;
147 exit(-1) ;
148 }
149 if (nt > n2c) {
150 cout << "cftcossins: nt > n2c : nt = " << nt << " , n2c = "
151 << n2c << endl ;
152 abort () ;
153 exit(-1) ;
154 }
155 if (n1f > n1c) {
156 cout << "cftcossins: n1f > n1c : n1f = " << n1f << " , n1c = "
157 << n1c << endl ;
158 abort () ;
159 exit(-1) ;
160 }
161 if (n3f > n3c) {
162 cout << "cftcossins: n3f > n3c : n3f = " << n3f << " , n3c = "
163 << n3c << endl ;
164 abort () ;
165 exit(-1) ;
166 }
167
168// Nombre de points pour la FFT:
169 int nm1 = nt - 1;
170 int nm1s2 = nm1 / 2;
171
172// Recherche des tables pour la FFT:
173 Tbl* pg = 0x0 ;
174 fftw_plan p = prepare_fft(nm1, pg) ;
175 Tbl& g = *pg ;
176
177// Recherche de la table des sin(psi) :
178 double* sinp = cheb_ini(nt);
179
180// boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
181// et 0 a dimf[2])
182
183 int n2n3f = n2f * n3f ;
184 int n2n3c = n2c * n3c ;
185
186//=======================================================================
187// Cas m pair
188//=======================================================================
189
190 j = 0 ;
191
192 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
193 // (car nul)
194
195//--------------------------------------------------------------------
196// partie cos(m phi) avec m pair : transformation en sin(l) theta)
197//--------------------------------------------------------------------
198
199 for (k=0; k<n3f; k++) {
200
201 int i0 = n2n3f * j + k ; // indice de depart
202 double* ff0 = ff + i0 ; // tableau des donnees a transformer
203
204 i0 = n2n3c * j + k ; // indice de depart
205 double* cf0 = cf + i0 ; // tableau resultat
206
207// Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
208//---------------------------------------------
209
210 for ( i = 1; i < nm1s2 ; i++ ) {
211 int isym = nm1 - i ;
212 int ix = n3f * i ;
213 int ixsym = n3f * isym ;
214// ... F+(theta)
215 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
216// ... F_(theta) sin(theta)
217 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
218 g.set(i) = fp + fms ;
219 g.set(isym) = fp - fms ;
220 }
221//... cas particuliers:
222 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
223 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
224
225// Developpement de G en series de Fourier par une FFT
226//----------------------------------------------------
227
228 fftw_execute(p) ;
229
230// Coefficients pairs du developmt. sin(l theta) de f
231//----------------------------------------------------
232// Ces coefficients sont egaux aux coefficients en sinus du developpement
233// de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
234// de fftw) :
235
236 double fac = -2. / double(nm1) ;
237 cf0[0] = 0. ;
238 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
239 cf0[n3c*nm1] = 0. ;
240
241// Coefficients impairs du developmt. en sin(l theta) de f
242//---------------------------------------------------------
243// 1. Coef. c'_k (recurrence amorcee a partir de zero):
244// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
245// (si fftw donnait reellement les coef. en sinus, il faudrait le
246// remplacer par un -2.)
247
248 cf0[n3c] = -fac * g(0);
249 fac *= -2. ;
250 for ( i = 3; i < nt; i += 2 ) {
251 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
252 }
253
254 } // fin de la boucle sur r
255
256//------------------------------------------------------------------------
257// partie sin(m phi) avec m pair : transformation en sin(l theta)
258//------------------------------------------------------------------------
259
260 j++ ;
261
262 if ( j != n1f-1 ) {
263// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
264// pas nuls
265
266 for (k=0; k<n3f; k++) {
267
268 int i0 = n2n3f * j + k ; // indice de depart
269 double* ff0 = ff + i0 ; // tableau des donnees a transformer
270
271 i0 = n2n3c * j + k ; // indice de depart
272 double* cf0 = cf + i0 ; // tableau resultat
273
274// Fonction G(theta) = F+(theta)sin(theta) + F_(theta)
275//---------------------------------------------
276 for ( i = 1; i < nm1s2 ; i++ ) {
277 int isym = nm1 - i ;
278 int ix = n3f * i ;
279 int ixsym = n3f * isym ;
280// ... F+(theta)
281 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
282// ... F_(theta) sin(theta)
283 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
284 g.set(i) = fp + fms ;
285 g.set(isym) = fp - fms ;
286 }
287//... cas particuliers:
288 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
289 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
290
291// Developpement de G en series de Fourier par une FFT
292//----------------------------------------------------
293
294 fftw_execute(p) ;
295
296// Coefficients pairs du developmt. sin(l theta) de f
297//----------------------------------------------------
298// Ces coefficients sont egaux aux coefficients en sinus du developpement
299// de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
300// de fftw) :
301
302 double fac = -2. / double(nm1) ;
303 cf0[0] = 0. ;
304 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(nm1 - i/2) ;
305 cf0[n3c*nm1] = 0. ;
306
307// Coefficients impairs du developmt. en sin(l theta) de f
308//---------------------------------------------------------
309// 1. Coef. c'_k (recurrence amorcee a partir de zero):
310// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
311// (si fftw donnait reellement les coef. en sinus, il faudrait le
312// remplacer par un -2.)
313
314 cf0[n3c] = -fac * g(0);
315 fac *= -2. ;
316 for ( i = 3; i < nt; i += 2 ) {
317 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(i/2) ;
318 }
319
320 } // fin de la boucle sur r
321
322 } // fin du cas ou le calcul etait necessaire (i.e. ou les
323 // coef en phi n'etaient pas nuls)
324
325
326// On passe au cas m impair suivant:
327 j+=3 ;
328
329 } // fin de la boucle sur les cas m pair
330
331 if (n1f<=3) // cas m=0 seulement (symetrie axiale)
332 return ;
333
334//=======================================================================
335// Cas m impair
336//=======================================================================
337
338 j = 2 ;
339
340 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
341 // (car nul)
342
343//------------------------------------------------------------------------
344// partie cos(m phi) avec m impair : transformation en cos(l theta)
345//------------------------------------------------------------------------
346
347
348 for (k=0; k<n3f; k++) {
349
350 int i0 = n2n3f * j + k ; // indice de depart
351 double* ff0 = ff + i0 ; // tableau des donnees a transformer
352
353 i0 = n2n3c * j + k ; // indice de depart
354 double* cf0 = cf + i0 ; // tableau resultat
355
356
357// Valeur en theta=0 de la partie antisymetrique de F, F_ :
358 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
359
360// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
361//---------------------------------------------
362 for ( i = 1; i < nm1s2 ; i++ ) {
363 int isym = nm1 - i ;
364 int ix = n3f * i ;
365 int ixsym = n3f * isym ;
366// ... F+(theta)
367 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
368// ... F_(theta) sin(psi)
369 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
370 g.set(i) = fp + fms ;
371 g.set(isym) = fp - fms ;
372 }
373//... cas particuliers:
374 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
375 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
376
377// Developpement de G en series de Fourier par une FFT
378//----------------------------------------------------
379
380 fftw_execute(p) ;
381
382// Coefficients pairs du developmt. cos(l theta) de f
383//----------------------------------------------------
384// Ces coefficients sont egaux aux coefficients en cosinus du developpement
385// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
386// de fftw) :
387
388 double fac = 2./double(nm1) ;
389 cf0[0] = g(0) / double(nm1) ;
390 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
391 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
392
393// Coefficients impairs du developmt. en cos(l theta) de f
394//---------------------------------------------------------
395// 1. Coef. c'_k (recurrence amorcee a partir de zero):
396// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
397// (si fftw donnait reellement les coef. en sinus, il faudrait le
398// remplacer par un -2.)
399 fac *= 2. ;
400 cf0[n3c] = 0 ;
401 double som = 0;
402 for ( i = 3; i < nt; i += 2 ) {
403 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
404 som += cf0[n3c*i] ;
405 }
406
407// 2. Calcul de c_1 :
408 double c1 = ( fmoins0 - som ) / nm1s2 ;
409
410// 3. Coef. c_k avec k impair:
411 cf0[n3c] = c1 ;
412 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
413
414
415 } // fin de la boucle sur r
416
417//--------------------------------------------------------------------
418// partie sin(m phi) avec m impair : transformation en cos(l theta)
419//--------------------------------------------------------------------
420
421 j++ ;
422
423 if ( (j != 1) && (j != n1f-1 ) ) {
424// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
425// pas nuls
426
427 for (k=0; k<n3f; k++) {
428
429 int i0 = n2n3f * j + k ; // indice de depart
430 double* ff0 = ff + i0 ; // tableau des donnees a transformer
431
432 i0 = n2n3c * j + k ; // indice de depart
433 double* cf0 = cf + i0 ; // tableau resultat
434
435
436// Valeur en theta=0 de la partie antisymetrique de F, F_ :
437 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
438
439// Fonction G(theta) = F+(theta) + F_(theta) sin(theta)
440//---------------------------------------------
441 for ( i = 1; i < nm1s2 ; i++ ) {
442 int isym = nm1 - i ;
443 int ix = n3f * i ;
444 int ixsym = n3f * isym ;
445// ... F+(theta)
446 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
447// ... F_(theta) sin(psi)
448 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
449 g.set(i) = fp + fms ;
450 g.set(isym) = fp - fms ;
451 }
452//... cas particuliers:
453 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
454 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
455
456// Developpement de G en series de Fourier par une FFT
457//----------------------------------------------------
458
459 fftw_execute(p) ;
460
461// Coefficients pairs du developmt. cos(l theta) de f
462//----------------------------------------------------
463// Ces coefficients sont egaux aux coefficients en cosinus du developpement
464// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
465// de fftw) :
466
467 double fac = 2./double(nm1) ;
468 cf0[0] = g(0) / double(nm1) ;
469 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
470 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
471
472// Coefficients impairs du developmt. en cos(l theta) de f
473//---------------------------------------------------------
474// 1. Coef. c'_k (recurrence amorcee a partir de zero):
475// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
476// (si fftw donnait reellement les coef. en sinus, il faudrait le
477// remplacer par un -2.)
478 fac *= 2. ;
479 cf0[n3c] = 0 ;
480 double som = 0;
481 for ( i = 3; i < nt; i += 2 ) {
482 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
483 som += cf0[n3c*i] ;
484 }
485
486// 2. Calcul de c_1 :
487 double c1 = ( fmoins0 - som ) / nm1s2 ;
488
489// 3. Coef. c_k avec k impair:
490 cf0[n3c] = c1 ;
491 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
492
493
494 } // fin de la boucle sur r
495 } // fin du cas ou le calcul etait necessaire (i.e. ou les
496 // coef en phi n'etaient pas nuls)
497
498// On passe au cas m impair suivant:
499 j+=3 ;
500
501 } // fin de la boucle sur les cas m pair
502
503}
504}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord sinp
Definition map.h:735