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