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