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