LORENE
FFTW3/cftcossinsp.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 en sin(2*l*theta) ou cos((2*l+1)*theta) (suivant la parite
27  * de l'indice m en phi) sur le deuxieme indice (theta)
28  * d'un tableau 3-D representant une fonction antisymetrique par rapport
29  * au plan z=0.
30  * Utilise la bibliotheque fftw.
31  *
32  * Entree:
33  * -------
34  * int* deg : tableau du nombre effectif de degres de liberte dans chacune
35  * des 3 dimensions: le nombre de points de collocation
36  * en theta est nt = deg[1] et doit etre de la forme
37  * nt = 2*p + 1
38  * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
39  * dimensions.
40  * On doit avoir dimf[1] >= deg[1] = nt.
41  *
42  * double* ff : tableau des valeurs de la fonction aux nt points de
43  * de collocation
44  *
45  * theta_l = pi/2 l/(nt-1) 0 <= l <= nt-1
46  *
47  * L'espace memoire correspondant a ce
48  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
49  * etre alloue avant l'appel a la routine.
50  * Les valeurs de la fonction doivent etre stokees
51  * dans le tableau ff comme suit
52  * f( theta_l ) = ff[ dimf[1]*dimf[2] * m + k + dimf[2] * l ]
53  * ou m et k sont les indices correspondant a
54  * phi et r respectivement.
55  * NB: cette routine suppose que la transformation en phi a deja ete
56  * effectuee: ainsi m est un indice de Fourier, non un indice de
57  * point de collocation en phi.
58  *
59  * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
60  * dimensions.
61  * On doit avoir dimc[1] >= deg[1] = nt.
62  * Sortie:
63  * -------
64  * double* cf : tableau des coefficients c_l de la fonction definis
65  * comme suit (a r et phi fixes)
66  *
67  * pour m pair:
68  * f(theta) = som_{l=0}^{nt-2} c_l sin( 2 l theta ) .
69  * pour m impair:
70  * f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) .
71  *
72  * L'espace memoire correspondant a ce
73  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
74  * etre alloue avant l'appel a la routine.
75  * Le coefficient c_l (0 <= l <= nt-1) est stoke dans
76  * le tableau cf comme suit
77  * c_l = cf[ dimc[1]*dimc[2] * m + k + dimc[2] * l ]
78  * ou m et k sont les indices correspondant a
79  * phi et r respectivement.
80  * Pour m pair, c_0 = c_{nt-1} = 0.
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: cftcossinsp.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
90  * $Log: cftcossinsp.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.0 1999/02/22 15:47:12 hyc
121  * *** empty log message ***
122  *
123  *
124  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinsp.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
125  *
126  */
127 
128 // headers du C
129 #include <cstdlib>
130 #include <fftw3.h>
131 
132 //Lorene prototypes
133 #include "tbl.h"
134 
135 // Prototypage des sous-routines utilisees:
136 namespace Lorene {
137 fftw_plan prepare_fft(int, Tbl*&) ;
138 double* cheb_ini(const int) ;
139 double* chebimp_ini(const int ) ;
140 //*****************************************************************************
141 
142 void cftcossinsp(const int* deg, const int* dimf, double* ff, const int* dimc,
143  double* cf)
144 {
145 
146 int i, j, k ;
147 
148 // Dimensions des tableaux ff et cf :
149  int n1f = dimf[0] ;
150  int n2f = dimf[1] ;
151  int n3f = dimf[2] ;
152  int n1c = dimc[0] ;
153  int n2c = dimc[1] ;
154  int n3c = dimc[2] ;
155 
156 // Nombre de degres de liberte en theta :
157  int nt = deg[1] ;
158 
159 // Tests de dimension:
160  if (nt > n2f) {
161  cout << "cftcossinsp: nt > n2f : nt = " << nt << " , n2f = "
162  << n2f << endl ;
163  abort () ;
164  exit(-1) ;
165  }
166  if (nt > n2c) {
167  cout << "cftcossinsp: nt > n2c : nt = " << nt << " , n2c = "
168  << n2c << endl ;
169  abort () ;
170  exit(-1) ;
171  }
172  if (n1f > n1c) {
173  cout << "cftcossinsp: n1f > n1c : n1f = " << n1f << " , n1c = "
174  << n1c << endl ;
175  abort () ;
176  exit(-1) ;
177  }
178  if (n3f > n3c) {
179  cout << "cftcossinsp: n3f > n3c : n3f = " << n3f << " , n3c = "
180  << n3c << endl ;
181  abort () ;
182  exit(-1) ;
183  }
184 
185 // Nombre de points pour la FFT:
186  int nm1 = nt - 1;
187  int nm1s2 = nm1 / 2;
188 
189 // Recherche des tables pour la FFT:
190  Tbl* pg = 0x0 ;
191  fftw_plan p = prepare_fft(nm1, pg) ;
192  Tbl& g = *pg ;
193 
194 // Recherche de la table des sin(psi) :
195  double* sinp = cheb_ini(nt);
196 
197 // Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}) :
198  double* x = chebimp_ini(nt);
199 
200 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
201 // et 0 a dimf[2])
202 
203  int n2n3f = n2f * n3f ;
204  int n2n3c = n2c * n3c ;
205 
206 //=======================================================================
207 // Cas m pair
208 //=======================================================================
209 
210  j = 0 ;
211 
212  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
213  // (car nul)
214 
215 //------------------------------------------------------------------------
216 // partie cos(m phi) avec m pair : transformation en sin((2 l) theta)
217 //------------------------------------------------------------------------
218 
219  for (k=0; k<n3f; k++) {
220 
221  int i0 = n2n3f * j + k ; // indice de depart
222  double* ff0 = ff + i0 ; // tableau des donnees a transformer
223 
224  i0 = n2n3c * j + k ; // indice de depart
225  double* cf0 = cf + i0 ; // tableau resultat
226 
227 /*
228  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
229  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
230  */
231 
232 // Fonction G(psi) = F+(psi) sin(psi) + F_(psi)
233 //---------------------------------------------
234  for ( i = 1; i < nm1s2 ; i++ ) {
235 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
236  int isym = nm1 - i ;
237 // ... indice (dans le tableau ff0) du point theta correspondant a psi
238  int ix = n3f * i ;
239 // ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
240  int ixsym = n3f * isym ;
241 // ... F+(psi) sin(psi)
242  double fps = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
243 // ... F_(psi)
244  double fm = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
245  g.set(i) = fps + fm ;
246  g.set(isym) = fps - fm ;
247  }
248 //... cas particuliers:
249  g.set(0) = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ) ;
250  g.set(nm1s2) = ff0[ n3f*nm1s2 ] ;
251 
252 // Developpement de G en series de Fourier par une FFT
253 //----------------------------------------------------
254  fftw_execute(p) ;
255 
256 // Coefficients pairs du developmt. sin(2l theta) de f
257 //----------------------------------------------------
258 // Ces coefficients sont egaux aux coefficients en sinus du developpement
259 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
260 // de fftw: si fftw donnait reellement les coefficients en sinus,
261 // il faudrait le remplacer par un +1) :
262  double fac = 2./double(nm1) ;
263  cf0[0] = 0. ;
264  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = -fac * g(nm1 - i/2);
265  cf0[n3c*nm1] = 0. ;
266 
267 // Coefficients impairs du developmt. en sin(2l theta) de f
268 //---------------------------------------------------------
269 // Ces coefficients sont obtenus par une recurrence a partir des
270 // coefficients en cosinus du developpement de G en series de Fourier
271 // (le facteur 4/nm1 vient de la normalisation
272 // de fftw: si fftw donnait reellement les coefficients en cosinus,
273 // il faudrait le remplacer par un +2.)
274 
275  cf0[n3c] = fac * g(0) ;
276  fac *= 2. ;
277  for ( i = 3; i < nt; i += 2 ) {
278  cf0[ n3c*i ] = cf0[ n3c*(i-2) ] + fac * g(i/2) ;
279  }
280 
281  } // fin de la boucle sur r
282 
283 //--------------------------------------------------------------------
284 // partie sin(m phi) avec m pair : transformation en sin((2 l) theta)
285 //--------------------------------------------------------------------
286 
287  j++ ;
288 
289  if ( (j != 1) && (j != n1f-1 ) ) {
290 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
291 // pas nuls
292 
293  for (k=0; k<n3f; k++) {
294 
295  int i0 = n2n3f * j + k ; // indice de depart
296  double* ff0 = ff + i0 ; // tableau des donnees a transformer
297 
298  i0 = n2n3c * j + k ; // indice de depart
299  double* cf0 = cf + i0 ; // tableau resultat
300 
301 /*
302  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
303  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
304  */
305 
306 // Fonction G(psi) = F+(psi) sin(psi) + F_(psi)
307 //---------------------------------------------
308  for ( i = 1; i < nm1s2 ; i++ ) {
309 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
310  int isym = nm1 - i ;
311 // ... indice (dans le tableau ff0) du point theta correspondant a psi
312  int ix = n3f * i ;
313 // ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
314  int ixsym = n3f * isym ;
315 // ... F+(psi) sin(psi)
316  double fps = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
317 // ... F_(psi)
318  double fm = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
319  g.set(i) = fps + fm ;
320  g.set(isym) = fps - fm ;
321  }
322 //... cas particuliers:
323  g.set(0) = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ) ;
324  g.set(nm1s2) = ff0[ n3f*nm1s2 ] ;
325 
326 // Developpement de G en series de Fourier par une FFT
327 //----------------------------------------------------
328  fftw_execute(p) ;
329 
330 // Coefficients pairs du developmt. sin(2l theta) de f
331 //----------------------------------------------------
332 // Ces coefficients sont egaux aux coefficients en sinus du developpement
333 // de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
334 // de fftw: si fftw donnait reellement les coefficients en sinus,
335 // il faudrait le remplacer par un +1) :
336  double fac = 2./double(nm1) ;
337  cf0[0] = 0. ;
338  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = -fac * g(nm1 - i/2);
339  cf0[n3c*nm1] = 0. ;
340 
341 // Coefficients impairs du developmt. en sin(2l theta) de f
342 //---------------------------------------------------------
343 // Ces coefficients sont obtenus par une recurrence a partir des
344 // coefficients en cosinus du developpement de G en series de Fourier
345 // (le facteur 4/nm1 vient de la normalisation
346 // de fftw: si fftw donnait reellement les coefficients en cosinus,
347 // il faudrait le remplacer par un +2.)
348 
349  cf0[n3c] = fac * g(0) ;
350  fac *= 2. ;
351  for ( i = 3; i < nt; i += 2 ) {
352  cf0[ n3c*i ] = cf0[ n3c*(i-2) ] + fac * g(i/2) ;
353  }
354 
355  } // fin de la boucle sur r
356 
357  } // fin du cas ou le calcul etait necessaire (i.e. ou les
358  // coef en phi n'etaient pas nuls)
359 
360 // On passe au cas m pair suivant:
361  j+=3 ;
362 
363  } // fin de la boucle sur les cas m pair
364 
365  if (n1f<=3) // cas m=0 seulement (symetrie axiale)
366  return ;
367 
368 //=======================================================================
369 // Cas m impair
370 //=======================================================================
371 
372  j = 2 ;
373 
374  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
375  // (car nul)
376 
377 //--------------------------------------------------------------------
378 // partie cos(m phi) avec m impair : transformation en cos((2 l+1) theta)
379 //--------------------------------------------------------------------
380 
381  for (k=0; k<n3f; k++) {
382 
383  int i0 = n2n3f * j + k ; // indice de depart
384  double* ff0 = ff + i0 ; // tableau des donnees a transformer
385 
386  i0 = n2n3c * j + k ; // indice de depart
387  double* cf0 = cf + i0 ; // tableau resultat
388 
389 // Multiplication de la fonction par x=cos(theta) (pour la rendre paire)
390 // NB: dans les commentaires qui suivent, on note h(x) = x f(x).
391 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
392 // tableau cf0).
393  for (i=0; i<nt-1; i++) cf0[n3c*i] = x[nm1-i] * ff0[n3f*i] ;
394  cf0[n3c*nm1] = 0 ;
395 
396 /*
397  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
398  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
399  */
400 
401 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
402  double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
403 
404 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
405 //---------------------------------------------
406  for ( i = 1; i < nm1s2 ; i++ ) {
407 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
408  int isym = nm1 - i ;
409 // ... indice (dans le tableau cf0) du point theta correspondant a psi
410  int ix = n3c * i ;
411 // ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
412  int ixsym = n3c * isym ;
413 // ... F+(psi)
414  double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
415 // ... F_(psi) sin(psi)
416  double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
417  g.set(i) = fp + fms ;
418  g.set(isym) = fp - fms ;
419  }
420 //... cas particuliers:
421  g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
422  g.set(nm1s2) = cf0[ n3c*nm1s2 ];
423 
424 // Developpement de G en series de Fourier par une FFT
425 //----------------------------------------------------
426 
427  fftw_execute(p) ;
428 
429 // Coefficients pairs du developmt. de Tchebyshev de h
430 //----------------------------------------------------
431 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
432 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
433 // de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
434 // remplacer par un +1.) :
435 
436  double fac = 2./double(nm1) ;
437  cf0[0] = g(0) / double(nm1) ;
438  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;
439  cf0[n3c*nm1] = g(nm1s2) ;
440 
441 // Coefficients impairs du developmt. de Tchebyshev de h
442 //------------------------------------------------------
443 // 1. Coef. c'_k (recurrence amorcee a partir de zero)
444 // Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
445 // (si fftw donnait reellement les coef. en sinus, il faudrait le
446 // remplacer par un -2.)
447  fac *= 2. ;
448  cf0[n3c] = 0 ;
449  double som = 0;
450  for ( i = 3; i < nt; i += 2 ) {
451  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
452  som += cf0[n3c*i] ;
453  }
454 
455 // 2. Calcul de c_1 :
456  double c1 = ( fmoins0 - som ) / nm1s2 ;
457 
458 // 3. Coef. c_k avec k impair:
459  cf0[n3c] = c1 ;
460  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
461 
462 
463 // Coefficients de f en fonction de ceux de h
464 //-------------------------------------------
465 
466  cf0[0] = 2* cf0[0] ;
467  for (i=1; i<nm1; i++) {
468  cf0[n3c*i] = 2 * cf0[n3c*i] - cf0[n3c*(i-1)] ;
469  }
470  cf0[n3c*nm1] = 0 ;
471 
472  } // fin de la boucle sur r
473 
474 //------------------------------------------------------------------------
475 // partie sin(m phi) avec m impair : transformation en cos((2 l+1) theta)
476 //------------------------------------------------------------------------
477 
478  j++ ;
479 
480  if ( j != n1f-1 ) {
481 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
482 // pas nuls
483 
484  for (k=0; k<n3f; k++) {
485 
486  int i0 = n2n3f * j + k ; // indice de depart
487  double* ff0 = ff + i0 ; // tableau des donnees a transformer
488 
489  i0 = n2n3c * j + k ; // indice de depart
490  double* cf0 = cf + i0 ; // tableau resultat
491 
492 // Multiplication de la fonction par x=cos(theta) (pour la rendre paire)
493 // NB: dans les commentaires qui suivent, on note h(x) = x f(x).
494 // (Les valeurs de h dans l'espace des configurations sont stokees dans le
495 // tableau cf0).
496  for (i=0; i<nt-1; i++) cf0[n3c*i] = x[nm1-i] * ff0[n3f*i] ;
497  cf0[n3c*nm1] = 0 ;
498 
499 /*
500  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
501  * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
502  */
503 
504 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
505  double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
506 
507 // Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
508 //---------------------------------------------
509  for ( i = 1; i < nm1s2 ; i++ ) {
510 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
511  int isym = nm1 - i ;
512 // ... indice (dans le tableau cf0) du point theta correspondant a psi
513  int ix = n3c * i ;
514 // ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
515  int ixsym = n3c * isym ;
516 // ... F+(psi)
517  double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
518 // ... F_(psi) sin(psi)
519  double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
520  g.set(i) = fp + fms ;
521  g.set(isym) = fp - fms ;
522  }
523 //... cas particuliers:
524  g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
525  g.set(nm1s2) = cf0[ n3c*nm1s2 ];
526 
527 // Developpement de G en series de Fourier par une FFT
528 //----------------------------------------------------
529 
530  fftw_execute(p) ;
531 
532 // Coefficients pairs du developmt. de Tchebyshev de h
533 //----------------------------------------------------
534 // Ces coefficients sont egaux aux coefficients en cosinus du developpement
535 // de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
536 // de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
537 // remplacer par un +1.) :
538 
539  double fac = 2./double(nm1) ;
540  cf0[0] = g(0) / double(nm1) ;
541  for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;
542  cf0[n3c*nm1] = g(nm1s2) ;
543 
544 // Coefficients impairs du developmt. de Tchebyshev de h
545 //------------------------------------------------------
546 // 1. Coef. c'_k (recurrence amorcee a partir de zero)
547 // Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
548 // (si fftw donnait reellement les coef. en sinus, il faudrait le
549 // remplacer par un -2.)
550  fac *= 2. ;
551  cf0[n3c] = 0 ;
552  double som = 0;
553  for ( i = 3; i < nt; i += 2 ) {
554  cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
555  som += cf0[n3c*i] ;
556  }
557 
558 // 2. Calcul de c_1 :
559  double c1 = ( fmoins0 - som ) / nm1s2 ;
560 
561 // 3. Coef. c_k avec k impair:
562  cf0[n3c] = c1 ;
563  for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
564 
565 
566 // Coefficients de f en fonction de ceux de h
567 //-------------------------------------------
568 
569  cf0[0] = 2* cf0[0] ;
570  for (i=1; i<nm1; i++) {
571  cf0[n3c*i] = 2 * cf0[n3c*i] - cf0[n3c*(i-1)] ;
572  }
573  cf0[n3c*nm1] = 0 ;
574 
575  } // fin de la boucle sur r
576 
577  } // fin du cas ou le calcul etait necessaire (i.e. ou les
578  // coef en phi n'etaient pas nuls)
579 
580 
581 // On passe au cas m impair suivant:
582  j+=3 ;
583 
584  } // fin de la boucle sur les cas m impair
585 
586 }
587 }
Lorene prototypes.
Definition: app_hor.h:67