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