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