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