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