LORENE
FFTW3/cfrchebpimi.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 de Tchebyshev T_{2k+1}/T_{2k} sur le troisieme indice (indice
27 * correspondant a r) d'un tableau 3-D decrivant par exemple d/dr d'une
28 * fonction symetrique par rapport au plan equatorial z = 0 et sans aucune
29 * autre symetrie, cad que l'on effectue
30 * 1/ un developpement en polynomes de Tchebyshev impairs pour m pair
31 * 2/ un developpement en polynomes de Tchebyshev pairs pour m impair
32 *
33 * Utilise la bibliotheque fftw.
34 *
35 *
36 * Entree:
37 * -------
38 * int* deg : tableau du nombre effectif de degres de liberte dans chacune
39 * des 3 dimensions: le nombre de points de collocation
40 * en r est nr = deg[2] et doit etre de la forme
41 * nr = 2*p + 1
42 * int* dimf : tableau du nombre d'elements de ff dans chacune des trois
43 * dimensions.
44 * On doit avoir dimf[2] >= deg[2] = nr.
45 *
46 * double* ff : tableau des valeurs de la fonction aux nr points de
47 * de collocation
48 *
49 * x_i = sin( pi/2 i/(nr-1) ) 0 <= i <= nr-1
50 *
51 * Les valeurs de la fonction doivent etre stokees dans le
52 * tableau ff comme suit
53 * f( x_i ) = ff[ dimf[1]*dimf[2] * j + dimf[2] * k + i ]
54 * ou j et k sont les indices correspondant a phi et theta
55 * respectivement.
56 * L'espace memoire correspondant a ce pointeur doit etre
57 * dimf[0]*dimf[1]*dimf[2] et doit etre alloue avant l'appel a
58 * la routine.
59 *
60 * int* dimc : tableau du nombre d'elements de cf dans chacune des trois
61 * dimensions.
62 * On doit avoir dimc[2] >= deg[2] = nr.
63 *
64 * Sortie:
65 * -------
66 * double* cf : tableau des nr coefficients c_i de la fonction definis
67 * comme suit (a theta et phi fixes)
68 *
69 * -- pour m pair (i.e. j = 0, 1, 4, 5, 8, 9, ...) :
70 *
71 * f(x) = som_{i=0}^{nr-2} c_i T_{2i+1}(x)
72 *
73 * ou T_{2i+1}(x) designe le polynome de Tchebyshev de
74 * degre 2i+1.
75 *
76 * -- pour m impair (i.e. j = 2, 3, 6, 7, 10, 11, ...) :
77 *
78 * f(x) = som_{i=0}^{nr-1} c_i T_{2i}(x) ,
79 *
80 * ou T_{2i}(x) designe le polynome de Tchebyshev de
81 * degre 2i.
82 *
83 * Les coefficients c_i (0 <= i <= nr-1) sont stokes dans
84 * le tableau cf comme suit
85 * c_i = cf[ dimc[1]*dimc[2] * j + dimc[2] * k + i ]
86 * ou j et k sont les indices correspondant a phi et theta
87 * respectivement.
88 * L'espace memoire correspondant a ce pointeur doit etre
89 * dimc[0]*dimc[1]*dimc[2] et doit avoir ete alloue avant
90 * l'appel a la routine.
91 *
92 * NB: Si le pointeur ff est egal a cf, la routine ne travaille que sur un
93 * seul tableau, qui constitue une entree/sortie.
94 */
95
96/*
97 * $Id: cfrchebpimi.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
98 * $Log: cfrchebpimi.C,v $
99 * Revision 1.4 2016/12/05 16:18:04 j_novak
100 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
101 *
102 * Revision 1.3 2014/10/13 08:53:18 j_novak
103 * Lorene classes and functions now belong to the namespace Lorene.
104 *
105 * Revision 1.2 2014/10/06 15:18:47 j_novak
106 * Modified #include directives to use c++ syntax.
107 *
108 * Revision 1.1 2004/12/21 17:06:02 j_novak
109 * Added all files for using fftw3.
110 *
111 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
112 * Suppressed the directive #include <malloc.h> for malloc is defined
113 * in <stdlib.h>
114 *
115 * Revision 1.3 2002/10/16 14:36:44 j_novak
116 * Reorganization of #include instructions of standard C++, in order to
117 * use experimental version 3 of gcc.
118 *
119 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
120 * Modification of declaration of Fortran 77 prototypes for
121 * a better portability (in particular on IBM AIX systems):
122 * All Fortran subroutine names are now written F77_* and are
123 * defined in the new file C++/Include/proto_f77.h.
124 *
125 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
126 * LORENE
127 *
128 * Revision 2.0 1999/02/22 15:48:21 hyc
129 * *** empty log message ***
130 *
131 *
132 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfrchebpimi.C,v 1.4 2016/12/05 16:18:04 j_novak Exp $
133 *
134 */
135
136// headers du C
137#include <cassert>
138#include <cstdlib>
139#include <fftw3.h>
140
141//Lorene prototypes
142#include "tbl.h"
143
144// Prototypage des sous-routines utilisees:
145namespace Lorene {
146fftw_plan prepare_fft(int, Tbl*&) ;
147double* cheb_ini(const int) ;
148double* chebimp_ini(const int ) ;
149
150//*****************************************************************************
151
152void cfrchebpimi(const int* deg, const int* dimf, double* ff, const int* dimc,
153 double* cf)
154
155{
156
157int i, j, k ;
158
159// Dimensions des tableaux ff et cf :
160 int n1f = dimf[0] ;
161 int n2f = dimf[1] ;
162 int n3f = dimf[2] ;
163 int n1c = dimc[0] ;
164 int n2c = dimc[1] ;
165 int n3c = dimc[2] ;
166
167// Nombres de degres de liberte en r :
168 int nr = deg[2] ;
169
170// Tests de dimension:
171 if (nr > n3f) {
172 cout << "cfrchebpimi: nr > n3f : nr = " << nr << " , n3f = "
173 << n3f << endl ;
174 abort () ;
175 exit(-1) ;
176 }
177 if (nr > n3c) {
178 cout << "cfrchebpimi: nr > n3c : nr = " << nr << " , n3c = "
179 << n3c << endl ;
180 abort () ;
181 exit(-1) ;
182 }
183 if (n1f > n1c) {
184 cout << "cfrchebpimi: n1f > n1c : n1f = " << n1f << " , n1c = "
185 << n1c << endl ;
186 abort () ;
187 exit(-1) ;
188 }
189 if (n2f > n2c) {
190 cout << "cfrchebpimi: n2f > n2c : n2f = " << n2f << " , n2c = "
191 << n2c << endl ;
192 abort () ;
193 exit(-1) ;
194 }
195
196// Nombre de points pour la FFT:
197 int nm1 = nr - 1;
198 int nm1s2 = nm1 / 2;
199
200// Recherche des tables pour la FFT:
201 Tbl* pg = 0x0 ;
202 fftw_plan p = prepare_fft(nm1, pg) ;
203 Tbl& g = *pg ;
204
205// Recherche de la table des sin(psi) :
206 double* sinp = cheb_ini(nr);
207
208// Recherche de la table des points de collocations x_k :
209 double* x = chebimp_ini(nr);
210
211// boucle sur phi et theta
212
213 int n2n3f = n2f * n3f ;
214 int n2n3c = n2c * n3c ;
215
216//=======================================================================
217// Cas m pair
218//=======================================================================
219
220 j = 0 ;
221
222 while (j<n1f-1) {
223
224//------------------------------------------------------------------------
225// partie cos(m phi) avec m pair : developpement en T_{2i+1}(x)
226//------------------------------------------------------------------------
227
228 for (k=0; k<n2f; k++) {
229
230 int i0 = n2n3f * j + n3f * k ; // indice de depart
231 double* ff0 = ff + i0 ; // tableau des donnees a transformer
232
233 i0 = n2n3c * j + n3c * k ; // indice de depart
234 double* cf0 = cf + i0 ; // tableau resultat
235
236// Multiplication de la fonction par x (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 cf0[0] = 0 ;
241 for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
242
243/*
244 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
245 * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)).
246 */
247
248// Valeur en psi=0 de la partie antisymetrique de F, F_ :
249 double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
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 x correspondant a psi
257 int ix = nm1 - i ;
258// ... indice (dans le tableau cf0) du point x correspondant a sym(psi)
259 int ixsym = nm1 - isym ;
260
261// ... F+(psi)
262 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
263// ... F_(psi) sin(psi)
264 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
265 g.set(i) = fp + fms ;
266 g.set(isym) = fp - fms ;
267 }
268//... cas particuliers:
269 g.set(0) = 0.5 * ( cf0[nm1] + cf0[0] );
270 g.set(nm1s2) = cf0[nm1s2];
271
272// Developpement de G en series de Fourier par une FFT
273//----------------------------------------------------
274
275 fftw_execute(p) ;
276
277// Coefficients pairs du developmt. de Tchebyshev de h
278//----------------------------------------------------
279// Ces coefficients sont egaux aux coefficients en cosinus du developpement
280// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
281// de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
282// remplacer par un +1.) :
283
284 double fac = 2./double(nm1) ;
285 cf0[0] = g(0) / double(nm1) ;
286 for (i=2; i<nm1; i += 2 ) cf0[i] = fac * g(i/2) ;
287 cf0[nm1] = g(nm1s2) / double(nm1) ;
288
289// Coefficients impairs du developmt. de Tchebyshev de h
290//------------------------------------------------------
291// 1. Coef. c'_k (recurrence amorcee a partir de zero)
292// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
293// (si fftw donnait reellement les coef. en sinus, il faudrait le
294// remplacer par un -2.)
295 fac *= 2 ;
296 cf0[1] = 0 ;
297 double som = 0;
298 for ( i = 3; i < nr; i += 2 ) {
299 cf0[i] = cf0[i-2] + fac * g(nm1 - i/2) ;
300 som += cf0[i] ;
301 }
302
303// 2. Calcul de c_1 :
304 double c1 = ( fmoins0 - som ) / nm1s2 ;
305
306// 3. Coef. c_k avec k impair:
307 cf0[1] = c1 ;
308 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
309
310// Coefficients de f en fonction de ceux de h
311//-------------------------------------------
312
313 cf0[0] = 2* cf0[0] ;
314 for (i=1; i<nm1; i++) {
315 cf0[i] = 2 * cf0[i] - cf0[i-1] ;
316 }
317 cf0[nm1] = 0 ;
318
319
320 } // fin de la boucle sur theta
321
322
323//------------------------------------------------------------------------
324// partie sin(m phi) avec m pair : developpement en T_{2i+1}(x)
325//------------------------------------------------------------------------
326
327 j++ ;
328
329 if ( (j != 1) && (j != n1f-1) ) {
330// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
331// pas nuls
332
333 for (k=0; k<n2f; k++) {
334
335 int i0 = n2n3f * j + n3f * k ; // indice de depart
336 double* ff0 = ff + i0 ; // tableau des donnees a transformer
337
338 i0 = n2n3c * j + n3c * k ; // indice de depart
339 double* cf0 = cf + i0 ; // tableau resultat
340
341// Multiplication de la fonction par x (pour la rendre paire)
342// NB: dans les commentaires qui suivent, on note h(x) = x f(x).
343// (Les valeurs de h dans l'espace des configurations sont stokees dans le
344// tableau cf0).
345 cf0[0] = 0 ;
346 for (i=1; i<nr; i++) cf0[i] = x[i] * ff0[i] ;
347
348/*
349 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
350 * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)).
351 */
352
353// Valeur en psi=0 de la partie antisymetrique de F, F_ :
354 double fmoins0 = 0.5 * ( cf0[nm1] - cf0[0] );
355
356// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
357//---------------------------------------------
358 for ( i = 1; i < nm1s2 ; i++ ) {
359// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
360 int isym = nm1 - i ;
361// ... indice (dans le tableau cf0) du point x correspondant a psi
362 int ix = nm1 - i ;
363// ... indice (dans le tableau cf0) du point x correspondant a sym(psi)
364 int ixsym = nm1 - isym ;
365
366// ... F+(psi)
367 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
368// ... F_(psi) sin(psi)
369 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
370 g.set(i) = fp + fms ;
371 g.set(isym) = fp - fms ;
372 }
373//... cas particuliers:
374 g.set(0) = 0.5 * ( cf0[nm1] + cf0[0] );
375 g.set(nm1s2) = cf0[nm1s2];
376
377// Developpement de G en series de Fourier par une FFT
378//----------------------------------------------------
379
380 fftw_execute(p) ;
381
382// Coefficients pairs du developmt. de Tchebyshev de h
383//----------------------------------------------------
384// Ces coefficients sont egaux aux coefficients en cosinus du developpement
385// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
386// de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
387// remplacer par un +1.) :
388
389 double fac = 2./double(nm1) ;
390 cf0[0] = g(0) / double(nm1) ;
391 for (i=2; i<nm1; i += 2 ) cf0[i] = fac * g(i/2) ;
392 cf0[nm1] = g(nm1s2) / double(nm1) ;
393
394// Coefficients impairs du developmt. de Tchebyshev de h
395//------------------------------------------------------
396// 1. Coef. c'_k (recurrence amorcee a partir de zero)
397// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
398// (si fftw donnait reellement les coef. en sinus, il faudrait le
399// remplacer par un -2.)
400 fac *= 2 ;
401 cf0[1] = 0 ;
402 double som = 0;
403 for ( i = 3; i < nr; i += 2 ) {
404 cf0[i] = cf0[i-2] + fac * g(nm1 - i/2) ;
405 som += cf0[i] ;
406 }
407
408// 2. Calcul de c_1 :
409 double c1 = ( fmoins0 - som ) / nm1s2 ;
410
411// 3. Coef. c_k avec k impair:
412 cf0[1] = c1 ;
413 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
414
415// Coefficients de f en fonction de ceux de h
416//-------------------------------------------
417
418 cf0[0] = 2* cf0[0] ;
419 for (i=1; i<nm1; i++) {
420 cf0[i] = 2 * cf0[i] - cf0[i-1] ;
421 }
422 cf0[nm1] = 0 ;
423
424 } // fin de la boucle sur theta
425
426 } // fin du cas ou le calcul etait necessaire (i.e. ou les
427 // coef en phi n'etaient pas nuls)
428
429// On passe au cas m pair suivant:
430 j+=3 ;
431
432 } // fin de la boucle sur les cas m pair
433
434 if (n1f<=3) // cas m=0 seulement (symetrie axiale)
435 return ;
436
437
438//=======================================================================
439// Cas m impair
440//=======================================================================
441
442 j = 2 ;
443
444 while (j<n1f-1) {
445
446//--------------------------------------------------------------------
447// partie cos(m phi) avec m impair : developpement en T_{2i}(x)
448//--------------------------------------------------------------------
449
450 for (k=0; k<n2f; k++) {
451
452 int i0 = n2n3f * j + n3f * k ; // indice de depart
453 double* ff0 = ff + i0 ; // tableau des donnees a transformer
454
455 i0 = n2n3c * j + n3c * k ; // indice de depart
456 double* cf0 = cf + i0 ; // tableau resultat
457
458/*
459 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
460 * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)).
461 */
462
463// Valeur en psi=0 de la partie antisymetrique de F, F_ :
464 double fmoins0 = 0.5 * ( ff0[nm1] - ff0[0] );
465
466// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
467//---------------------------------------------
468 for ( i = 1; i < nm1s2 ; i++ ) {
469// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
470 int isym = nm1 - i ;
471// ... indice (dans le tableau ff0) du point x correspondant a psi
472 int ix = nm1 - i ;
473// ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
474 int ixsym = nm1 - isym ;
475
476// ... F+(psi)
477 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
478// ... F_(psi) sin(psi)
479 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
480 g.set(i) = fp + fms ;
481 g.set(isym) = fp - fms ;
482 }
483//... cas particuliers:
484 g.set(0) = 0.5 * ( ff0[nm1] + ff0[0] );
485 g.set(nm1s2) = ff0[nm1s2];
486
487// Developpement de G en series de Fourier par une FFT
488//----------------------------------------------------
489
490 fftw_execute(p) ;
491
492// Coefficients pairs du developmt. de Tchebyshev de f
493//----------------------------------------------------
494// Ces coefficients sont egaux aux coefficients en cosinus du developpement
495// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
496// de fftw) :
497
498 double fac = 2./double(nm1) ;
499 cf0[0] = g(0) / double(nm1) ;
500 for (i=2; i<nm1; i += 2) cf0[i] = fac*g(i/2) ;
501 cf0[nm1] = g(nm1s2) / double(nm1) ;
502
503// Coefficients impairs du developmt. de Tchebyshev de f
504//------------------------------------------------------
505// 1. Coef. c'_k (recurrence amorcee a partir de zero)
506// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
507// (si fftw donnait reellement les coef. en sinus, il faudrait le
508// remplacer par un -2.)
509 fac *= 2. ;
510 cf0[1] = 0. ;
511 double som = 0.;
512 for ( i = 3; i < nr; i += 2 ) {
513 cf0[i] = cf0[i-2] + fac * g(nm1-i/2) ;
514 som += cf0[i] ;
515 }
516
517// 2. Calcul de c_1 :
518 double c1 = ( fmoins0 - som ) / nm1s2 ;
519
520// 3. Coef. c_k avec k impair:
521 cf0[1] = c1 ;
522 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
523
524 } // fin de la boucle sur theta
525
526
527//--------------------------------------------------------------------
528// partie sin(m phi) avec m impair : developpement en T_{2i}(x)
529//--------------------------------------------------------------------
530
531 j++ ;
532
533 if ( j != n1f-1 ) {
534// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
535// pas nuls
536
537 for (k=0; k<n2f; k++) {
538
539 int i0 = n2n3f * j + n3f * k ; // indice de depart
540 double* ff0 = ff + i0 ; // tableau des donnees a transformer
541
542 i0 = n2n3c * j + n3c * k ; // indice de depart
543 double* cf0 = cf + i0 ; // tableau resultat
544
545/*
546 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
547 * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)).
548 */
549
550// Valeur en psi=0 de la partie antisymetrique de F, F_ :
551 double fmoins0 = 0.5 * ( ff0[nm1] - ff0[0] );
552
553// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
554//---------------------------------------------
555 for ( i = 1; i < nm1s2 ; i++ ) {
556// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
557 int isym = nm1 - i ;
558// ... indice (dans le tableau ff0) du point x correspondant a psi
559 int ix = nm1 - i ;
560// ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
561 int ixsym = nm1 - isym ;
562
563// ... F+(psi)
564 double fp = 0.5 * ( ff0[ix] + ff0[ixsym] ) ;
565// ... F_(psi) sin(psi)
566 double fms = 0.5 * ( ff0[ix] - ff0[ixsym] ) * sinp[i] ;
567 g.set(i) = fp + fms ;
568 g.set(isym) = fp - fms ;
569 }
570//... cas particuliers:
571 g.set(0) = 0.5 * ( ff0[nm1] + ff0[0] );
572 g.set(nm1s2) = ff0[nm1s2];
573
574// Developpement de G en series de Fourier par une FFT
575//----------------------------------------------------
576
577 fftw_execute(p) ;
578
579// Coefficients pairs du developmt. de Tchebyshev de f
580//----------------------------------------------------
581// Ces coefficients sont egaux aux coefficients en cosinus du developpement
582// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
583// de fftw) :
584
585 double fac = 2./double(nm1) ;
586 cf0[0] = g(0) / double(nm1) ;
587 for (i=2; i<nm1; i += 2) cf0[i] = fac*g(i/2) ;
588 cf0[nm1] = g(nm1s2) / double(nm1) ;
589
590// Coefficients impairs du developmt. de Tchebyshev de f
591//------------------------------------------------------
592// 1. Coef. c'_k (recurrence amorcee a partir de zero)
593// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
594// (si fftw donnait reellement les coef. en sinus, il faudrait le
595// remplacer par un -2.)
596 fac *= 2. ;
597 cf0[1] = 0. ;
598 double som = 0.;
599 for ( i = 3; i < nr; i += 2 ) {
600 cf0[i] = cf0[i-2] + fac * g(nm1-i/2) ;
601 som += cf0[i] ;
602 }
603
604// 2. Calcul de c_1 :
605 double c1 = ( fmoins0 - som ) / nm1s2 ;
606
607// 3. Coef. c_k avec k impair:
608 cf0[1] = c1 ;
609 for ( i = 3; i < nr; i += 2 ) cf0[i] += c1 ;
610
611 } // fin de la boucle sur theta
612
613 } // fin du cas ou le calcul etait necessaire (i.e. ou les
614 // coef en phi n'etaient pas nuls)
615
616// On passe au cas m impair suivant:
617 j+=3 ;
618
619 } // fin de la boucle sur les cas m impair
620
621}
622}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord x
x coordinate centered on the grid
Definition map.h:738
Coord sinp
Definition map.h:735