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