LORENE
FFTW3/citcossinc.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/*
27 * Transformation inverse cos(l*theta) ou sin(l*theta) (suivant la
28 * parite de l'indice m en phi) sur le deuxieme indice (theta)
29 * d'un tableau 3-D representant une fonction symetrique par rapport
30 * au plan z=0.
31 * Utilise la bibliotheque fftw.
32 *
33 * Entree:
34 * -------
35 * int* deg : tableau du nombre effectif de degres de liberte dans chacune
36 * des 3 dimensions: le nombre de points de collocation
37 * en theta est nt = deg[1] et doit etre de la forme
38 * nt = 2*p + 1
39 * int* dimc : tableau du nombre d'elements de cc dans chacune des trois
40 * dimensions.
41 * On doit avoir dimc[1] >= deg[1] = nt.
42 *
43 * double* cf : tableau des coefficients c_l de la fonction definis
44 * comme suit (a r et phi fixes)
45 *
46 * pour m pair:
47 * f(theta) = som_{l=0}^{nt-1} c_l cos( l theta ) .
48 * pour m impair:
49 * f(theta) = som_{l=0}^{nt-2} c_l sin( l theta ) .
50 *
51 * L'espace memoire correspondant a ce
52 * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
53 * avoir ete alloue avant l'appel a la routine.
54 * Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
55 * le tableau cf comme suit
56 * c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
57 * ou j et k sont les indices correspondant a
58 * phi et r respectivement.
59 *
60 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
61 * dimensions.
62 * On doit avoir dimf[1] >= deg[1] = nt.
63 *
64 * Sortie:
65 * -------
66 * double* ff : tableau des valeurs de la fonction aux nt points de
67 * de collocation
68 *
69 * theta_l = pi l/(nt-1) 0 <= l <= nt-1
70 *
71 * L'espace memoire correspondant a ce
72 * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
73 * avoir ete alloue avant l'appel a la routine.
74 * Les valeurs de la fonction sont stokees
75 * dans le tableau ff comme suit
76 * f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
77 * ou j et k sont les indices correspondant a
78 * phi et r respectivement.
79 *
80 * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un
81 * seul tableau, qui constitue une entree/sortie.
82 *
83 */
84
85/*
86 * $Id: citcossinc.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
87 * $Log: citcossinc.C,v $
88 * Revision 1.4 2016/12/05 16:18:05 j_novak
89 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
90 *
91 * Revision 1.3 2014/10/13 08:53:20 j_novak
92 * Lorene classes and functions now belong to the namespace Lorene.
93 *
94 * Revision 1.2 2014/10/06 15:18:50 j_novak
95 * Modified #include directives to use c++ syntax.
96 *
97 * Revision 1.1 2004/12/21 17:06:03 j_novak
98 * Added all files for using fftw3.
99 *
100 * Revision 1.1 2004/11/23 15:13:50 m_forot
101 * Added the bases for the cases without any equatorial symmetry
102 * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
103 *
104 *
105 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinc.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
106 *
107 */
108// headers du C
109#include <cstdlib>
110#include <fftw3.h>
111
112//Lorene prototypes
113#include "tbl.h"
114
115// Prototypage des sous-routines utilisees:
116namespace Lorene {
117fftw_plan back_fft(int, Tbl*&) ;
118double* cheb_ini(const int) ;
119//*****************************************************************************
120
121void citcossinc(const int* deg, const int* dimc, double* cf, const int* dimf,
122 double* ff)
123{
124
125int i, j, k ;
126
127// Dimensions des tableaux ff et cf :
128 int n1f = dimf[0] ;
129 int n2f = dimf[1] ;
130 int n3f = dimf[2] ;
131 int n1c = dimc[0] ;
132 int n2c = dimc[1] ;
133 int n3c = dimc[2] ;
134
135// Nombres de degres de liberte en theta :
136 int nt = deg[1] ;
137
138// Tests de dimension:
139 if (nt > n2f) {
140 cout << "citcossinc: nt > n2f : nt = " << nt << " , n2f = "
141 << n2f << endl ;
142 abort () ;
143 exit(-1) ;
144 }
145 if (nt > n2c) {
146 cout << "citcossinc: nt > n2c : nt = " << nt << " , n2c = "
147 << n2c << endl ;
148 abort () ;
149 exit(-1) ;
150 }
151 if (n1c > n1f) {
152 cout << "citcossinc: n1c > n1f : n1c = " << n1c << " , n1f = "
153 << n1f << endl ;
154 abort () ;
155 exit(-1) ;
156 }
157 if (n3c > n3f) {
158 cout << "citcossinc: n3c > n3f : n3c = " << n3c << " , n3f = "
159 << n3f << endl ;
160 abort () ;
161 exit(-1) ;
162 }
163
164// Nombre de points pour la FFT:
165 int nm1 = nt - 1;
166 int nm1s2 = nm1 / 2;
167
168// Recherche des tables pour la FFT:
169 Tbl* pg = 0x0 ;
170 fftw_plan p = back_fft(nm1, pg) ;
171 Tbl& g = *pg ;
172
173// Recherche de la table des sin(psi) :
174 double* sinp = cheb_ini(nt);
175
176// boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
177// et 0 a dimf[2])
178
179 int n2n3f = n2f * n3f ;
180 int n2n3c = n2c * n3c ;
181
182//=======================================================================
183// Cas m pair
184//=======================================================================
185
186 j = 0 ;
187
188 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
189 // (car nul)
190
191//-----------------------------------------------------------------------
192// partie cos(m phi) avec m pair : transformation cos( l theta) inverse
193//-----------------------------------------------------------------------
194
195 for (k=0; k<n3c; k++) {
196
197 int i0 = n2n3c * j + k ; // indice de depart
198 double* cf0 = cf + i0 ; // tableau des donnees a transformer
199
200 i0 = n2n3f * j + k ; // indice de depart
201 double* ff0 = ff + i0 ; // tableau resultat
202
203
204// Coefficients impairs de G
205//--------------------------
206
207 double c1 = cf0[n3c] ;
208
209 double som = 0;
210 ff0[n3f] = 0 ;
211 for ( i = 3; i < nt; i += 2 ) {
212 ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
213 som += ff0[ n3f*i ] ;
214 }
215
216// Valeur en theta=0 de la partie antisymetrique de F, F_ :
217 double fmoins0 = nm1s2 * c1 + som ;
218
219// Coef. impairs de G
220// NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
221// donnait exactement les coef. des sinus, ce facteur serait -0.5.
222 for ( i = 3; i < nt; i += 2 ) {
223 g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
224 }
225
226
227// Coefficients pairs de G
228//------------------------
229// Ces coefficients sont egaux aux coefficients pairs du developpement de
230// f.
231// NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
232// donnait exactement les coef. des cosinus, ce facteur serait 1.
233
234 g.set(0) = cf0[0] ;
235 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;
236 g.set(nm1s2) = cf0[ n3c*nm1 ] ;
237
238// Transformation de Fourier inverse de G
239//---------------------------------------
240
241// FFT inverse
242 fftw_execute(p) ;
243
244// Valeurs de f deduites de celles de G
245//-------------------------------------
246
247 for ( i = 1; i < nm1s2 ; i++ ) {
248// ... indice du pt symetrique de psi par rapport a pi/2:
249 int isym = nm1 - i ;
250
251 double fp = 0.5 * ( g(i) + g(isym) ) ;
252 double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
253 ff0[ n3f*i ] = fp + fm ;
254 ff0[ n3f*isym ] = fp - fm ;
255 }
256
257//... cas particuliers:
258 ff0[0] = g(0) + fmoins0 ;
259 ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
260 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
261 } // fin de la boucle sur r
262
263//-----------------------------------------------------------------------
264// partie sin(m phi) avec m pair : transformation cos(l theta) inverse
265//-----------------------------------------------------------------------
266
267 j++ ;
268
269 if ( (j != 1) && (j != n1f-1 ) ) {
270// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
271// pas nuls
272
273 for (k=0; k<n3c; k++) {
274
275 int i0 = n2n3c * j + k ; // indice de depart
276 double* cf0 = cf + i0 ; // tableau des donnees a transformer
277
278 i0 = n2n3f * j + k ; // indice de depart
279 double* ff0 = ff + i0 ; // tableau resultat
280
281// Coefficients impairs de G
282//--------------------------
283
284 double c1 = cf0[n3c] ;
285
286 double som = 0;
287 ff0[n3f] = 0 ;
288 for ( i = 3; i < nt; i += 2 ) {
289 ff0[ n3f*i ] = cf0[ n3c*i ] - c1 ;
290 som += ff0[ n3f*i ] ;
291 }
292
293// Valeur en theta=0 de la partie antisymetrique de F, F_ :
294 double fmoins0 = nm1s2 * c1 + som ;
295
296// Coef. impairs de G
297// NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
298// donnait exactement les coef. des sinus, ce facteur serait -0.5.
299 for ( i = 3; i < nt; i += 2 ) {
300 g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
301 }
302
303
304// Coefficients pairs de G
305//------------------------
306// Ces coefficients sont egaux aux coefficients pairs du developpement de
307// f.
308// NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
309// donnait exactement les coef. des cosinus, ce facteur serait 1.
310
311 g.set(0) = cf0[0] ;
312 for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * cf0[ n3c*2*i ] ;
313 g.set(nm1s2) = cf0[ n3c*nm1 ] ;
314
315// Transformation de Fourier inverse de G
316//---------------------------------------
317
318// FFT inverse
319 fftw_execute(p) ;
320
321// Valeurs de f deduites de celles de G
322//-------------------------------------
323
324 for ( i = 1; i < nm1s2 ; i++ ) {
325// ... indice du pt symetrique de psi par rapport a pi/2:
326 int isym = nm1 - i ;
327
328 double fp = 0.5 * ( g(i) + g(isym) ) ;
329 double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
330 ff0[ n3f*i ] = fp + fm ;
331 ff0[ n3f*isym ] = fp - fm ;
332 }
333
334//... cas particuliers:
335 ff0[0] = g(0) + fmoins0 ;
336 ff0[ n3f*nm1 ] = g(0) - fmoins0 ;
337 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
338 } // fin de la boucle sur r
339
340 } // fin du cas ou le calcul etait necessaire (i.e. ou les
341 // coef en phi n'etaient pas nuls)
342
343// On passe au cas m pair suivant:
344 j+=3 ;
345
346 } // fin de la boucle sur les cas m pair
347
348//##
349 if (n1f<=3) // cas m=0 seulement (symetrie axiale)
350 return ;
351
352//=======================================================================
353// Cas m impair
354//=======================================================================
355
356 j = 2 ;
357
358 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
359 // (car nul)
360
361//--------------------------------------------------------------------------
362// partie cos(m phi) avec m impair : transformation sin(l theta) inv.
363//--------------------------------------------------------------------------
364
365 for (k=0; k<n3c; k++) {
366
367 int i0 = n2n3c * j + k ; // indice de depart
368 double* cf0 = cf + i0 ; // tableau des donnees a transformer
369
370 i0 = n2n3f * j + k ; // indice de depart
371 double* ff0 = ff + i0 ; // tableau resultat
372
373// Coefficients impairs de G
374//--------------------------
375
376 for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
377
378// Coefficients pairs de G
379//------------------------
380
381 g.set(0) = .5 * cf0[n3c] ;
382 for ( i = 3; i < nt; i += 2 ) {
383 g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
384 }
385 g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
386
387// Transformation de Fourier inverse de G
388//---------------------------------------
389
390// FFT inverse
391 fftw_execute(p) ;
392
393// Valeurs de f deduites de celles de G
394//-------------------------------------
395
396 for ( i = 1; i < nm1s2 ; i++ ) {
397// ... indice du pt symetrique de psi par rapport a pi/2:
398 int isym = nm1 - i ;
399
400 double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
401 double fm = 0.5 * ( g(i) - g(isym) ) ;
402 ff0[ n3f*i ] = fp + fm ;
403 ff0[ n3f*isym ] = fp - fm ;
404 }
405
406//... cas particuliers:
407 ff0[0] = 0. ;
408 ff0[ n3f*nm1 ] = -2*g(0) ;
409 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
410 } // fin de la boucle sur r
411
412//--------------------------------------------------------------------------
413// partie sin(m phi) avec m impair : transformation sin(l theta) inv.
414//--------------------------------------------------------------------------
415
416 j++ ;
417
418 if ( j != n1f-1 ) {
419// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
420// pas nuls
421
422 for (k=0; k<n3c; k++) {
423
424 int i0 = n2n3c * j + k ; // indice de depart
425 double* cf0 = cf + i0 ; // tableau des donnees a transformer
426
427 i0 = n2n3f * j + k ; // indice de depart
428 double* ff0 = ff + i0 ; // tableau resultat
429
430// Coefficients impairs de G
431//--------------------------
432
433 for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = -0.5 * cf0[ n3c*i ] ;
434
435// Coefficients pairs de G
436//------------------------
437
438 g.set(0) = .5 * cf0[n3c] ;
439 for ( i = 3; i < nt; i += 2 ) {
440 g.set(i/2) = .25 * ( cf0[ n3c*i ] - cf0[ n3c*(i-2) ] ) ;
441 }
442 g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
443
444// Transformation de Fourier inverse de G
445//---------------------------------------
446
447// FFT inverse
448 fftw_execute(p) ;
449
450// Valeurs de f deduites de celles de G
451//-------------------------------------
452
453 for ( i = 1; i < nm1s2 ; i++ ) {
454// ... indice du pt symetrique de psi par rapport a pi/2:
455 int isym = nm1 - i ;
456
457 double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
458 double fm = 0.5 * ( g(i) - g(isym) ) ;
459 ff0[ n3f*i ] = fp + fm ;
460 ff0[ n3f*isym ] = fp - fm ;
461 }
462
463//... cas particuliers:
464 ff0[0] = 0. ;
465 ff0[ n3f*nm1 ] = -2*g(0) ;
466 ff0[ n3f*nm1s2 ] = g(nm1s2) ;
467 } // fin de la boucle sur r
468
469 } // fin du cas ou le calcul etait necessaire (i.e. ou les
470 // coef en phi n'etaient pas nuls)
471
472// On passe au cas m impair suivant:
473 j+=3 ;
474
475 } // fin de la boucle sur les cas m impair
476
477}
478}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord sinp
Definition map.h:735