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