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