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