LORENE
FFTW3/citcossinci.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 * Transformation inverse cos((2*l+1)*theta) ou sin(2*l*theta) (suivant la
27 * parite de l'indice m en phi) sur le deuxieme indice (theta)
28 * d'un tableau 3-D representant une fonction antisymetrique 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 cf dans chacune des trois
39 * dimensions.
40 * On doit avoir dimc[1] >= deg[1] = nt.
41 *
42 * double* cf : tableau des coefficients c_l de la fonction definis
43 * comme suit (a r et phi fixes)
44 *
45 * pour m pair:
46 * f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) .
47 * pour m impair:
48 * f(theta) = som_{l=0}^{nt-2} c_l sin( 2 l theta ) .
49 *
50 * L'espace memoire correspondant a ce
51 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
52 * avoir ete alloue avant l'appel a la routine.
53 * Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
54 * le tableau cf comme suit
55 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
56 * ou j et k sont les indices correspondant a
57 * phi et r respectivement.
58 *
59 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
60 * dimensions.
61 * On doit avoir dimf[1] >= deg[1] = nt.
62 *
63 * Sortie:
64 * -------
65 * double* ff : tableau des valeurs de la fonction aux nt points de
66 * de collocation
67 *
68 * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
69 *
70 * L'espace memoire correspondant a ce
71 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
72 * avoir ete alloue avant l'appel a la routine.
73 * Les valeurs de la fonction sont stokees
74 * dans le tableau ff comme suit
75 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
76 * ou j et k sont les indices correspondant a
77 * phi et r respectivement.
78 *
79 * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
80 * seul tableau, qui constitue une entree/sortie.
81 *
82 */
83
84/*
85 * $Id: citcossinci.C,v 1.6 2016/12/05 16:18:05 j_novak Exp $
86 * $Log: citcossinci.C,v $
87 * Revision 1.6 2016/12/05 16:18:05 j_novak
88 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
89 *
90 * Revision 1.5 2014/10/13 08:53:20 j_novak
91 * Lorene classes and functions now belong to the namespace Lorene.
92 *
93 * Revision 1.4 2014/10/06 15:18:50 j_novak
94 * Modified #include directives to use c++ syntax.
95 *
96 * Revision 1.3 2013/04/25 15:46:06 j_novak
97 * Added special treatment in the case np = 1, for type_p = NONSYM.
98 *
99 * Revision 1.2 2004/12/27 14:27:28 j_novak
100 * Added forgotten "delete [] t1"
101 *
102 * Revision 1.1 2004/12/21 17:06:03 j_novak
103 * Added all files for using fftw3.
104 *
105 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
106 * Suppressed the directive #include <malloc.h> for malloc is defined
107 * in <stdlib.h>
108 *
109 * Revision 1.3 2002/10/16 14:36:53 j_novak
110 * Reorganization of #include instructions of standard C++, in order to
111 * use experimental version 3 of gcc.
112 *
113 * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
114 * Modification of declaration of Fortran 77 prototypes for
115 * a better portability (in particular on IBM AIX systems):
116 * All Fortran subroutine names are now written F77_* and are
117 * defined in the new file C++/Include/proto_f77.h.
118 *
119 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
120 * LORENE
121 *
122 * Revision 2.0 1999/02/22 15:42:35 hyc
123 * *** empty log message ***
124 *
125 *
126 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinci.C,v 1.6 2016/12/05 16:18:05 j_novak Exp $
127 *
128 */
129
130// headers du C
131#include <cstdlib>
132#include <fftw3.h>
133
134//Lorene prototypes
135#include "tbl.h"
136
137// Prototypage des sous-routines utilisees:
138namespace Lorene {
139fftw_plan back_fft(int, Tbl*&) ;
140double* cheb_ini(const int) ;
141double* chebimp_ini(const int ) ;
142//*****************************************************************************
143
144void citcossinci(const int* deg, const int* dimc, double* cf, const int* dimf,
145 double* ff)
146{
147
148int i, j, k ;
149
150// Dimensions des tableaux ff et cf :
151 int n1f = dimf[0] ;
152 int n2f = dimf[1] ;
153 int n3f = dimf[2] ;
154 int n1c = dimc[0] ;
155 int n2c = dimc[1] ;
156 int n3c = dimc[2] ;
157
158// Nombres de degres de liberte en theta :
159 int nt = deg[1] ;
160
161// Tests de dimension:
162 if (nt > n2f) {
163 cout << "citcossinci: nt > n2f : nt = " << nt << " , n2f = "
164 << n2f << endl ;
165 abort () ;
166 }
167 if (nt > n2c) {
168 cout << "citcossinci: nt > n2c : nt = " << nt << " , n2c = "
169 << n2c << endl ;
170 abort () ;
171 }
172 if ( (n1f > 1) && (n1c > n1f) ) {
173 cout << "citcossinci: n1c > n1f : n1c = " << n1c << " , n1f = "
174 << n1f << endl ;
175 abort () ;
176 }
177 if (n3c > n3f) {
178 cout << "citcossinci: n3c > n3f : n3c = " << n3c << " , n3f = "
179 << n3f << endl ;
180 abort () ;
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 double* t1 = new double[nt] ;
192
193// Recherche de la table des sin(psi) :
194 double* sinp = cheb_ini(nt);
195
196// Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}):
197 double* x = chebimp_ini(nt) ;
198
199// boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
200// et 0 a dimf[2])
201
202 int n2n3f = n2f * n3f ;
203 int n2n3c = n2c * n3c ;
204 int borne_phi = n1f-1 ;
205 if (n1f == 1) borne_phi = 1 ;
206
207//=======================================================================
208// Cas m pair
209//=======================================================================
210
211 j = 0 ;
212
213 while (j < borne_phi) { //le dernier coef en phi n'est pas considere
214 // (car nul)
215
216//-----------------------------------------------------------------------
217// partie cos(m phi) avec m pair : transformation cos((2 l+1) theta) inverse
218//-----------------------------------------------------------------------
219
220 for (k=0; k<n3c; k++) {
221
222 int i0 = n2n3c * j + k ; // indice de depart
223 double* cf0 = cf + i0 ; // tableau des donnees a transformer
224
225 i0 = n2n3f * j + k ; // indice de depart
226 double* ff0 = ff + i0 ; // tableau resultat
227
228// Calcul des coefficients du developpement en T_{2i}(x) de la fonction
229// h(x) := x f(x) (x=cos(theta)) a partir des coefficients de f
230// (resultat stoke dans le tableau t1 :
231 t1[0] = .5 * cf0[0] ;
232 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
233 t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
234
235/*
236 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
237 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
238 */
239
240// Calcul des coefficients de Fourier de la fonction
241// G(psi) = F+(psi) + F_(psi) sin(psi)
242// en fonction des coefficients en cos(2l theta) de f:
243
244// Coefficients impairs de G
245//--------------------------
246
247 double c1 = t1[1] ;
248
249 double som = 0;
250 ff0[n3f] = 0 ;
251 for ( i = 3; i < nt; i += 2 ) {
252 ff0[ n3f*i ] = t1[i] - c1 ;
253 som += ff0[ n3f*i ] ;
254 }
255
256// Valeur en psi=0 de la partie antisymetrique de F, F_ :
257 double fmoins0 = nm1s2 * c1 + som ;
258
259// Coef. impairs de G
260// NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
261// donnait exactement les coef. des sinus, ce facteur serait -0.5.
262 for ( i = 3; i < nt; i += 2 ) {
263 g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
264 }
265
266
267// Coefficients pairs de G
268//------------------------
269// Ces coefficients sont egaux aux coefficients pairs du developpement de
270// f.
271// NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
272// donnait exactement les coef. des cosinus, ce facteur serait 1.
273
274 g.set(0) = t1[0] ;
275 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
276 g.set(nm1s2) = t1[nm1] ;
277
278// Transformation de Fourier inverse de G
279//---------------------------------------
280
281// FFT inverse
282 fftw_execute(p) ;
283
284// Valeurs de f deduites de celles de G
285//-------------------------------------
286
287 for ( i = 1; i < nm1s2 ; i++ ) {
288// ... indice du pt symetrique de psi par rapport a pi/2:
289 int isym = nm1 - i ;
290
291 double fp = 0.5 * ( g(i) + g(isym) ) ;
292 double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
293 ff0[ n3f*i ] = ( fp + fm ) / x[isym] ;
294 ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
295 }
296
297//... cas particuliers:
298 ff0[0] = g(0) + fmoins0 ;
299 ff0[ n3f*nm1 ] = 0 ;
300 ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
301
302 } // fin de la boucle sur r
303
304//-----------------------------------------------------------------------
305// partie sin(m phi) avec m pair : transformation cos((2 l+1) theta) inverse
306//-----------------------------------------------------------------------
307
308 j++ ;
309
310 if ( (j != 1) && (j != borne_phi ) ) {
311// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
312// pas nuls
313
314 for (k=0; k<n3c; k++) {
315
316 int i0 = n2n3c * j + k ; // indice de depart
317 double* cf0 = cf + i0 ; // tableau des donnees a transformer
318
319 i0 = n2n3f * j + k ; // indice de depart
320 double* ff0 = ff + i0 ; // tableau resultat
321
322// Calcul des coefficients du developpement en T_{2i}(x) de la fonction
323// h(x) := x f(x) (x=cos(theta)) a partir des coefficients de f
324// (resultat stoke dans le tableau t1 :
325 t1[0] = .5 * cf0[0] ;
326 for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
327 t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
328
329/*
330 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
331 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
332 */
333
334// Calcul des coefficients de Fourier de la fonction
335// G(psi) = F+(psi) + F_(psi) sin(psi)
336// en fonction des coefficients en cos(2l theta) de f:
337
338// Coefficients impairs de G
339//--------------------------
340
341 double c1 = t1[1] ;
342
343 double som = 0;
344 ff0[n3f] = 0 ;
345 for ( i = 3; i < nt; i += 2 ) {
346 ff0[ n3f*i ] = t1[i] - c1 ;
347 som += ff0[ n3f*i ] ;
348 }
349
350// Valeur en psi=0 de la partie antisymetrique de F, F_ :
351 double fmoins0 = nm1s2 * c1 + som ;
352
353// Coef. impairs de G
354// NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
355// donnait exactement les coef. des sinus, ce facteur serait -0.5.
356 for ( i = 3; i < nt; i += 2 ) {
357 g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
358 }
359
360
361// Coefficients pairs de G
362//------------------------
363// Ces coefficients sont egaux aux coefficients pairs du developpement de
364// f.
365// NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
366// donnait exactement les coef. des cosinus, ce facteur serait 1.
367
368 g.set(0) = t1[0] ;
369 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ;
370 g.set(nm1s2) = t1[nm1] ;
371
372// Transformation de Fourier inverse de G
373//---------------------------------------
374
375// FFT inverse
376 fftw_execute(p) ;
377
378// Valeurs de f deduites de celles de G
379//-------------------------------------
380
381 for ( i = 1; i < nm1s2 ; i++ ) {
382// ... indice du pt symetrique de psi par rapport a pi/2:
383 int isym = nm1 - i ;
384
385 double fp = 0.5 * ( g(i) + g(isym) ) ;
386 double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
387 ff0[ n3f*i ] = ( fp + fm ) / x[isym] ;
388 ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
389 }
390
391//... cas particuliers:
392 ff0[0] = g(0) + fmoins0 ;
393 ff0[ n3f*nm1 ] = 0 ;
394 ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
395
396
397 } // fin de la boucle sur r
398
399 } // fin du cas ou le calcul etait necessaire (i.e. ou les
400 // coef en phi n'etaient pas nuls)
401
402// On passe au cas m pair suivant:
403 j+=3 ;
404
405 } // fin de la boucle sur les cas m pair
406
407//##
408 if (n1f<=3) {
409 delete [] t1 ;
410 return ;
411 }
412
413//=======================================================================
414// Cas m impair
415//=======================================================================
416
417 j = 2 ;
418
419 while (j < borne_phi) { //le dernier coef en phi n'est pas considere
420 // (car nul)
421
422//--------------------------------------------------------------------------
423// partie cos(m phi) avec m impair : transformation sin((2 l) theta) inv.
424//--------------------------------------------------------------------------
425
426 for (k=0; k<n3c; k++) {
427
428 int i0 = n2n3c * j + k ; // indice de depart
429 double* cf0 = cf + i0 ; // tableau des donnees a transformer
430
431 i0 = n2n3f * j + k ; // indice de depart
432 double* ff0 = ff + i0 ; // tableau resultat
433
434
435/*
436 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
437 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
438 */
439
440
441// Calcul des coefficients de Fourier de la fonction
442// G(psi) = F+(psi) + F_(psi) sin(psi)
443// en fonction des coefficients en sin(2l theta) de f:
444
445// Coefficients en sinus de G
446//---------------------------
447// Ces coefficients sont egaux aux coefficients pairs du developmt. en
448// sin(2l theta) de f (le facteur -.5 vient de la normalisation
449// de fftw: si fftw donnait reellement les coefficients en sinus,
450// il faudrait le remplacer par un +1) :
451 for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;
452
453// Coefficients en cosinus de G
454//-----------------------------
455// Ces coefficients se deduisent des coefficients impairs du developmt.
456// en sin(2l theta) de f (le facteur +.25 vient de la normalisation
457// de fftw: si fftw donnait reellement les coefficients en cosinus,
458// il faudrait le remplacer par un +.5)
459
460 g.set(0) = .5 * cf0[n3c] ;
461 for ( i = 1; i < nm1s2; i++ ) {
462 g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
463 }
464 g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
465
466
467// Transformation de Fourier inverse de G
468//---------------------------------------
469
470// FFT inverse
471 fftw_execute(p) ;
472
473// Valeurs de f deduites de celles de G
474//-------------------------------------
475
476 for ( i = 1; i < nm1s2 ; i++ ) {
477// ... indice du pt symetrique de psi par rapport a pi/2:
478 int isym = nm1 - i ;
479
480 double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
481 double fm = 0.5 * ( g(i) - g(isym) ) ;
482 ff0[ n3f*i ] = fp + fm ;
483 ff0[ n3f*isym ] = fp - fm ;
484 }
485
486//... cas particuliers:
487 ff0[0] = 0. ;
488 ff0[ n3f*nm1 ] = -2. * g(0) ;
489 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
490
491 } // fin de la boucle sur r
492
493
494//--------------------------------------------------------------------------
495// partie sin(m phi) avec m impair : transformation sin( (2 l) theta) inv.
496//--------------------------------------------------------------------------
497
498 j++ ;
499
500 if ( j != borne_phi ) {
501// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
502// pas nuls
503
504 for (k=0; k<n3c; k++) {
505
506 int i0 = n2n3c * j + k ; // indice de depart
507 double* cf0 = cf + i0 ; // tableau des donnees a transformer
508
509 i0 = n2n3f * j + k ; // indice de depart
510 double* ff0 = ff + i0 ; // tableau resultat
511
512
513/*
514 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
515 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
516 */
517
518
519// Calcul des coefficients de Fourier de la fonction
520// G(psi) = F+(psi) + F_(psi) sin(psi)
521// en fonction des coefficients en sin(2l theta) de f:
522
523// Coefficients en sinus de G
524//---------------------------
525// Ces coefficients sont egaux aux coefficients pairs du developmt. en
526// sin(2l theta) de f (le facteur -.5 vient de la normalisation
527// de fftw: si fftw donnait reellement les coefficients en sinus,
528// il faudrait le remplacer par un +1) :
529 for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;
530
531// Coefficients en cosinus de G
532//-----------------------------
533// Ces coefficients se deduisent des coefficients impairs du developmt.
534// en sin(2l theta) de f (le facteur +.25 vient de la normalisation
535// de fftw: si fftw donnait reellement les coefficients en cosinus,
536// il faudrait le remplacer par un +.5)
537
538 g.set(0) = .5 * cf0[n3c] ;
539 for ( i = 1; i < nm1s2; i++ ) {
540 g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
541 }
542 g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
543
544
545// Transformation de Fourier inverse de G
546//---------------------------------------
547
548// FFT inverse
549 fftw_execute(p) ;
550
551// Valeurs de f deduites de celles de G
552//-------------------------------------
553
554 for ( i = 1; i < nm1s2 ; i++ ) {
555// ... indice du pt symetrique de psi par rapport a pi/2:
556 int isym = nm1 - i ;
557
558 double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
559 double fm = 0.5 * ( g(i) - g(isym) ) ;
560 ff0[ n3f*i ] = fp + fm ;
561 ff0[ n3f*isym ] = fp - fm ;
562 }
563
564//... cas particuliers:
565 ff0[0] = 0. ;
566 ff0[ n3f*nm1 ] = -2. * g(0) ;
567 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
568
569 } // fin de la boucle sur r
570
571
572 } // fin du cas ou le calcul etait necessaire (i.e. ou les
573 // coef en phi n'etaient pas nuls)
574
575// On passe au cas m impair suivant:
576 j+=3 ;
577
578 } // fin de la boucle sur les cas m impair
579
580 delete [] t1 ;
581
582}
583}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord x
x coordinate centered on the grid
Definition map.h:738
Coord sinp
Definition map.h:735