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