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