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