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