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