LORENE
FFTW3/cftcossinsi.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 cftcossinsi_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinsi.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $" ;
24
25/*
26 * Transformation en sin((2*l+1)*theta) ou cos(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 symetrique par rapport
29 * au plan z=0.
30 * Utilise la bibiotheque 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+1) theta ) .
69 * pour m impair:
70 * f(theta) = som_{l=0}^{nt-1} c_l cos( 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 *
82 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
83 * seul tableau, qui constitue une entree/sortie.
84 *
85 */
86
87/*
88 * $Id: cftcossinsi.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $
89 * $Log: cftcossinsi.C,v $
90 * Revision 1.3 2014/10/13 08:53:19 j_novak
91 * Lorene classes and functions now belong to the namespace Lorene.
92 *
93 * Revision 1.2 2014/10/06 15:18:48 j_novak
94 * Modified #include directives to use c++ syntax.
95 *
96 * Revision 1.1 2004/12/21 17:06:02 j_novak
97 * Added all files for using fftw3.
98 *
99 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
100 * Suppressed the directive #include <malloc.h> for malloc is defined
101 * in <stdlib.h>
102 *
103 * Revision 1.3 2002/10/16 14:36:51 j_novak
104 * Reorganization of #include instructions of standard C++, in order to
105 * use experimental version 3 of gcc.
106 *
107 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
108 * Modification of declaration of Fortran 77 prototypes for
109 * a better portability (in particular on IBM AIX systems):
110 * All Fortran subroutine names are now written F77_* and are
111 * defined in the new file C++/Include/proto_f77.h.
112 *
113 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
114 * LORENE
115 *
116 * Revision 2.1 2000/01/27 12:16:13 eric
117 * Modif commentaires.
118 *
119 * Revision 2.0 1999/02/22 15:47:20 hyc
120 * *** empty log message ***
121 *
122 *
123 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinsi.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $
124 *
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:
136namespace Lorene {
137fftw_plan prepare_fft(int, Tbl*&) ;
138double* cheb_ini(const int) ;
139double* chebimp_ini(const int ) ;
140//*****************************************************************************
141
142void cftcossinsi(const int* deg, const int* dimf, double* ff, const int* dimc,
143 double* cf)
144{
145
146int 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 << "cftcossinsi: nt > n2f : nt = " << nt << " , n2f = "
162 << n2f << endl ;
163 abort () ;
164 exit(-1) ;
165 }
166 if (nt > n2c) {
167 cout << "cftcossinsi: nt > n2c : nt = " << nt << " , n2c = "
168 << n2c << endl ;
169 abort () ;
170 exit(-1) ;
171 }
172 if (n1f > n1c) {
173 cout << "cftcossinsi: n1f > n1c : n1f = " << n1f << " , n1c = "
174 << n1c << endl ;
175 abort () ;
176 exit(-1) ;
177 }
178 if (n3f > n3c) {
179 cout << "cftcossinsi: 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 sin( theta_l ) :
198 double* sinth = 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+1) 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// Multiplication de la fonction par sin(theta) (pour la rendre developpable
228// en cos(2l theta) )
229// NB: dans les commentaires qui suivent, on note
230// h(theta) = f(theta) sin(theta).
231// (Les valeurs de h dans l'espace des configurations sont stokees dans le
232// tableau cf0).
233 cf0[0] = 0 ;
234 for (i=1; i<nt; i++) cf0[n3c*i] = sinth[i] * ff0[n3f*i] ;
235
236/*
237 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
238 * reliee a theta par psi = 2 theta et F(psi) = h(theta(psi)).
239 */
240
241// Valeur en psi=0 de la partie antisymetrique de F, F_ :
242 double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
243
244// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
245//---------------------------------------------
246 for ( i = 1; i < nm1s2 ; i++ ) {
247// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
248 int isym = nm1 - i ;
249// ... indice (dans le tableau cf0) du point theta correspondant a psi
250 int ix = n3c * i ;
251// ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
252 int ixsym = n3c * isym ;
253// ... F+(psi)
254 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
255// ... F_(psi) sin(psi)
256 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
257 g.set(i) = fp + fms ;
258 g.set(isym) = fp - fms ;
259 }
260//... cas particuliers:
261 g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
262 g.set(nm1s2) = cf0[ n3c*nm1s2 ];
263
264// Developpement de G en series de Fourier par une FFT
265//----------------------------------------------------
266
267 fftw_execute(p) ;
268
269// Coefficients pairs du developmt. cos(2l theta) de h
270//----------------------------------------------------
271// Ces coefficients sont egaux aux coefficients en cosinus du developpement
272// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
273// de fftw) :
274
275 double fac = 2./double(nm1) ;
276 cf0[0] = g(0)/double(nm1) ;
277 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;
278 cf0[n3c*nm1] = g(nm1s2)/double(nm1) ;
279
280// Coefficients impairs du developmt. en cos(2l theta) de h
281//---------------------------------------------------------
282// 1. Coef. c'_k (recurrence amorcee a partir de zero):
283// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
284// (si fftw donnait reellement les coef. en sinus, il faudrait le
285// remplacer par un -2.)
286 fac *= 2. ;
287 cf0[n3c] = 0 ;
288 double som = 0;
289 for ( i = 3; i < nt; i += 2 ) {
290 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
291 som += cf0[n3c*i] ;
292 }
293
294// 2. Calcul de c_1 :
295 double c1 = ( fmoins0 - som ) / nm1s2 ;
296
297// 3. Coef. c_k avec k impair:
298 cf0[n3c] = c1 ;
299 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
300
301// Coefficients de f en fonction de ceux de h
302//-------------------------------------------
303
304 cf0[0] = 2* cf0[0] ;
305 for (i=1; i<nm1; i++) {
306 cf0[n3c*i] = 2 * cf0[n3c*i] + cf0[n3c*(i-1)] ;
307 }
308 cf0[n3c*nm1] = 0 ;
309
310 } // fin de la boucle sur r
311
312//--------------------------------------------------------------------
313// partie sin(m phi) avec m pair : transformation en sin((2*l+1) theta)
314//--------------------------------------------------------------------
315
316 j++ ;
317
318 if ( (j != 1) && (j != n1f-1 ) ) {
319// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
320// pas nuls
321
322 for (k=0; k<n3f; k++) {
323
324 int i0 = n2n3f * j + k ; // indice de depart
325 double* ff0 = ff + i0 ; // tableau des donnees a transformer
326
327 i0 = n2n3c * j + k ; // indice de depart
328 double* cf0 = cf + i0 ; // tableau resultat
329
330// Multiplication de la fonction par sin(theta) (pour la rendre developpable
331// en cos(2l theta) )
332// NB: dans les commentaires qui suivent, on note
333// h(theta) = f(theta) sin(theta).
334// (Les valeurs de h dans l'espace des configurations sont stokees dans le
335// tableau cf0).
336 cf0[0] = 0 ;
337 for (i=1; i<nt; i++) cf0[n3c*i] = sinth[i] * ff0[n3f*i] ;
338
339/*
340 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
341 * reliee a theta par psi = 2 theta et F(psi) = h(theta(psi)).
342 */
343
344// Valeur en psi=0 de la partie antisymetrique de F, F_ :
345 double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
346
347// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
348//---------------------------------------------
349 for ( i = 1; i < nm1s2 ; i++ ) {
350// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
351 int isym = nm1 - i ;
352// ... indice (dans le tableau cf0) du point theta correspondant a psi
353 int ix = n3c * i ;
354// ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
355 int ixsym = n3c * isym ;
356// ... F+(psi)
357 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
358// ... F_(psi) sin(psi)
359 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
360 g.set(i) = fp + fms ;
361 g.set(isym) = fp - fms ;
362 }
363//... cas particuliers:
364 g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
365 g.set(nm1s2) = cf0[ n3c*nm1s2 ];
366
367// Developpement de G en series de Fourier par une FFT
368//----------------------------------------------------
369
370 fftw_execute(p) ;
371
372// Coefficients pairs du developmt. cos(2l theta) de h
373//----------------------------------------------------
374// Ces coefficients sont egaux aux coefficients en cosinus du developpement
375// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
376// de fftw) :
377
378 double fac = 2./double(nm1) ;
379 cf0[0] = g(0)/double(nm1) ;
380 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;
381 cf0[n3c*nm1] = g(nm1s2)/double(nm1) ;
382
383// Coefficients impairs du developmt. en cos(2l theta) de h
384//---------------------------------------------------------
385// 1. Coef. c'_k (recurrence amorcee a partir de zero):
386// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
387// (si fftw donnait reellement les coef. en sinus, il faudrait le
388// remplacer par un -2.)
389 fac *= 2. ;
390 cf0[n3c] = 0 ;
391 double som = 0;
392 for ( i = 3; i < nt; i += 2 ) {
393 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
394 som += cf0[n3c*i] ;
395 }
396
397// 2. Calcul de c_1 :
398 double c1 = ( fmoins0 - som ) / nm1s2 ;
399
400// 3. Coef. c_k avec k impair:
401 cf0[n3c] = c1 ;
402 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
403
404// Coefficients de f en fonction de ceux de h
405//-------------------------------------------
406
407 cf0[0] = 2* cf0[0] ;
408 for (i=1; i<nm1; i++) {
409 cf0[n3c*i] = 2 * cf0[n3c*i] + cf0[n3c*(i-1)] ;
410 }
411 cf0[n3c*nm1] = 0 ;
412
413 } // fin de la boucle sur r
414
415 } // fin du cas ou le calcul etait necessaire (i.e. ou les
416 // coef en phi n'etaient pas nuls)
417
418// On passe au cas m pair suivant:
419 j+=3 ;
420
421 } // fin de la boucle sur les cas m pair
422
423 if (n1f<=3) // cas m=0 seulement (symetrie axiale)
424 return ;
425
426//=======================================================================
427// Cas m impair
428//=======================================================================
429
430 j = 2 ;
431
432 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
433 // (car nul)
434
435//--------------------------------------------------------------------
436// partie cos(m phi) avec m impair : transformation en cos(2 l theta)
437//--------------------------------------------------------------------
438
439 for (k=0; k<n3f; k++) {
440
441 int i0 = n2n3f * j + k ; // indice de depart
442 double* ff0 = ff + i0 ; // tableau des donnees a transformer
443
444 i0 = n2n3c * j + k ; // indice de depart
445 double* cf0 = cf + i0 ; // tableau resultat
446
447/*
448 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
449 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
450 */
451
452// Valeur en psi=0 de la partie antisymetrique de F, F_ :
453 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
454
455// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
456//---------------------------------------------
457 for ( i = 1; i < nm1s2 ; i++ ) {
458// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
459 int isym = nm1 - i ;
460// ... indice (dans le tableau ff0) du point theta correspondant a psi
461 int ix = n3f * i ;
462// ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
463 int ixsym = n3f * isym ;
464// ... F+(psi)
465 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
466// ... F_(psi) sin(psi)
467 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
468 g.set(i) = fp + fms ;
469 g.set(isym) = fp - fms ;
470 }
471//... cas particuliers:
472 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
473 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
474
475// Developpement de G en series de Fourier par une FFT
476//----------------------------------------------------
477
478 fftw_execute(p) ;
479
480// Coefficients pairs du developmt. cos(2l theta) de f
481//----------------------------------------------------
482// Ces coefficients sont egaux aux coefficients en cosinus du developpement
483// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
484// de fftw) :
485
486 double fac = 2./double(nm1) ;
487 cf0[0] = g(0) / double(nm1) ;
488 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
489 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
490
491// Coefficients impairs du developmt. en cos(2l theta) de f
492//---------------------------------------------------------
493// 1. Coef. c'_k (recurrence amorcee a partir de zero):
494// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
495// (si fftw donnait reellement les coef. en sinus, il faudrait le
496// remplacer par un -2.)
497 fac *= 2. ;
498 cf0[n3c] = 0 ;
499 double som = 0;
500 for ( i = 3; i < nt; i += 2 ) {
501 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
502 som += cf0[n3c*i] ;
503 }
504
505// 2. Calcul de c_1 :
506 double c1 = ( fmoins0 - som ) / nm1s2 ;
507
508// 3. Coef. c_k avec k impair:
509 cf0[n3c] = c1 ;
510 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
511
512
513 } // fin de la boucle sur r
514
515//------------------------------------------------------------------------
516// partie sin(m phi) avec m impair : transformation en cos(2 l theta)
517//------------------------------------------------------------------------
518
519 j++ ;
520
521 if ( j != n1f-1 ) {
522// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
523// pas nuls
524
525 for (k=0; k<n3f; k++) {
526
527 int i0 = n2n3f * j + k ; // indice de depart
528 double* ff0 = ff + i0 ; // tableau des donnees a transformer
529
530 i0 = n2n3c * j + k ; // indice de depart
531 double* cf0 = cf + i0 ; // tableau resultat
532
533/*
534 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
535 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
536 */
537
538// Valeur en psi=0 de la partie antisymetrique de F, F_ :
539 double fmoins0 = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] );
540
541// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
542//---------------------------------------------
543 for ( i = 1; i < nm1s2 ; i++ ) {
544// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
545 int isym = nm1 - i ;
546// ... indice (dans le tableau ff0) du point theta correspondant a psi
547 int ix = n3f * i ;
548// ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
549 int ixsym = n3f * isym ;
550// ... F+(psi)
551 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
552// ... F_(psi) sin(psi)
553 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
554 g.set(i) = fp + fms ;
555 g.set(isym) = fp - fms ;
556 }
557//... cas particuliers:
558 g.set(0) = 0.5 * ( ff0[0] + ff0[ n3f*nm1 ] );
559 g.set(nm1s2) = ff0[ n3f*nm1s2 ];
560
561// Developpement de G en series de Fourier par une FFT
562//----------------------------------------------------
563
564 fftw_execute(p) ;
565
566// Coefficients pairs du developmt. cos(2l theta) de f
567//----------------------------------------------------
568// Ces coefficients sont egaux aux coefficients en cosinus du developpement
569// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
570// de fftw) :
571
572 double fac = 2./double(nm1) ;
573 cf0[0] = g(0) / double(nm1) ;
574 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac*g(i/2) ;
575 cf0[n3c*nm1] = g(nm1s2) / double(nm1) ;
576
577// Coefficients impairs du developmt. en cos(2l theta) de f
578//---------------------------------------------------------
579// 1. Coef. c'_k (recurrence amorcee a partir de zero):
580// Le 4/nm1 en facteur de g[i] est du a la normalisation de fftw
581// (si fftw donnait reellement les coef. en sinus, il faudrait le
582// remplacer par un -2.)
583 fac *= 2. ;
584 cf0[n3c] = 0 ;
585 double som = 0;
586 for ( i = 3; i < nt; i += 2 ) {
587 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
588 som += cf0[n3c*i] ;
589 }
590
591// 2. Calcul de c_1 :
592 double c1 = ( fmoins0 - som ) / nm1s2 ;
593
594// 3. Coef. c_k avec k impair:
595 cf0[n3c] = c1 ;
596 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
597
598
599 } // fin de la boucle sur r
600
601 } // fin du cas ou le calcul etait necessaire (i.e. ou les
602 // coef en phi n'etaient pas nuls)
603
604
605// On passe au cas m impair suivant:
606 j+=3 ;
607
608 } // fin de la boucle sur les cas m impair
609
610}
611}
Lorene prototypes.
Definition app_hor.h:64