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