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