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 
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.5 2016/12/05 16:18:04 j_novak Exp $
100  * $Log: circhebpimp.C,v $
101  * Revision 1.5 2016/12/05 16:18:04 j_novak
102  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
103  *
104  * Revision 1.4 2014/10/15 12:48:21 j_novak
105  * Corrected namespace declaration.
106  *
107  * Revision 1.3 2014/10/13 08:53:17 j_novak
108  * Lorene classes and functions now belong to the namespace Lorene.
109  *
110  * Revision 1.2 2014/10/06 15:18:46 j_novak
111  * Modified #include directives to use c++ syntax.
112  *
113  * Revision 1.1 2004/12/21 17:06:01 j_novak
114  * Added all files for using fftw3.
115  *
116  * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
117  * Suppressed the directive #include <malloc.h> for malloc is defined
118  * in <stdlib.h>
119  *
120  * Revision 1.3 2002/10/16 14:36:53 j_novak
121  * Reorganization of #include instructions of standard C++, in order to
122  * use experimental version 3 of gcc.
123  *
124  * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
125  * Modification of declaration of Fortran 77 prototypes for
126  * a better portability (in particular on IBM AIX systems):
127  * All Fortran subroutine names are now written F77_* and are
128  * defined in the new file C++/Include/proto_f77.h.
129  *
130  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
131  * LORENE
132  *
133  * Revision 2.0 1999/02/22 15:43:10 hyc
134  * *** empty log message ***
135  *
136  *
137  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/circhebpimp.C,v 1.5 2016/12/05 16:18:04 j_novak Exp $
138  *
139  */
140 
141 
142 // headers du C
143 #include <cassert>
144 #include <cstdlib>
145 
146 // Prototypes of F77 subroutines
147 #include "headcpp.h"
148 #include "proto_f77.h"
149 
150 // Prototypage des sous-routines utilisees:
151 namespace Lorene {
152 int* facto_ini(int ) ;
153 double* trigo_ini(int ) ;
154 double* cheb_ini(const int) ;
155 double* chebimp_ini(const int ) ;
156 //*****************************************************************************
157 
158 void circhebpimp(const int* deg, const int* dimc, double* cf,
159  const int* dimf, double* ff)
160 
161 {
162 
163 int i, j, k ;
164 
165 // Dimensions des tableaux ff et cf :
166  int n1f = dimf[0] ;
167  int n2f = dimf[1] ;
168  int n3f = dimf[2] ;
169  int n1c = dimc[0] ;
170  int n2c = dimc[1] ;
171  int n3c = dimc[2] ;
172 
173 // Nombres de degres de liberte en r :
174  int nr = deg[2] ;
175 
176 // Tests de dimension:
177  if (nr > n3c) {
178  cout << "circhebpimp: nr > n3c : nr = " << nr << " , n3c = "
179  << n3c << endl ;
180  abort () ;
181  exit(-1) ;
182  }
183  if (nr > n3f) {
184  cout << "circhebpimp: nr > n3f : nr = " << nr << " , n3f = "
185  << n3f << endl ;
186  abort () ;
187  exit(-1) ;
188  }
189  if (n1c > n1f) {
190  cout << "circhebpimp: n1c > n1f : n1c = " << n1c << " , n1f = "
191  << n1f << endl ;
192  abort () ;
193  exit(-1) ;
194  }
195  if (n2c > n2f) {
196  cout << "circhebpimp: n2c > n2f : n2c = " << n2c << " , n2f = "
197  << n2f << endl ;
198  abort () ;
199  exit(-1) ;
200  }
201 
202 // Nombre de points pour la FFT:
203  int nm1 = nr - 1;
204  int nm1s2 = nm1 / 2;
205 
206 // Recherche des tables pour la FFT:
207  int* facto = facto_ini(nm1) ;
208  double* trigo = trigo_ini(nm1) ;
209 
210 // Recherche de la table des sin(psi) :
211  double* sinp = cheb_ini(nr);
212 
213 // Recherche de la table des points de collocations x_k :
214  double* x = chebimp_ini(nr);
215 
216  // tableau de travail t1 et g
217  // (la dimension nm1+2 = nr+1 est exigee par la routine fft991)
218  double* g = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
219  double* t1 = (double*)( malloc( (nm1+2)*sizeof(double) ) ) ;
220 
221 // Parametres pour la routine FFT991
222  int jump = 1 ;
223  int inc = 1 ;
224  int lot = 1 ;
225  int isign = 1 ;
226 
227 // boucle sur phi et theta
228 
229  int n2n3f = n2f * n3f ;
230  int n2n3c = n2c * n3c ;
231 
232 //=======================================================================
233 // Cas m pair
234 //=======================================================================
235 
236  j = 0 ;
237 
238  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
239  // (car nul)
240 
241 //--------------------------------------------------------------------
242 // partie cos(m phi) avec m pair : developpement en T_{2i}(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 /*
254  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
255  * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)).
256  */
257 
258 // Calcul des coefficients de Fourier de la fonction
259 // G(psi) = F+(psi) + F_(psi) sin(psi)
260 // en fonction des coefficients de Tchebyshev de f:
261 
262 // Coefficients impairs de G
263 //--------------------------
264 
265  double c1 = cf0[1] ;
266 
267  double som = 0;
268  ff0[1] = 0 ;
269  for ( i = 3; i < nr; i += 2 ) {
270  ff0[i] = cf0[i] - c1 ;
271  som += ff0[i] ;
272  }
273 
274 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
275  double fmoins0 = nm1s2 * c1 + som ;
276 
277 // Coef. impairs de G
278 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
279 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
280  g[1] = 0 ;
281  for ( i = 3; i < nr; i += 2 ) {
282  g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
283  }
284  g[nr] = 0 ;
285 
286 
287 // Coefficients pairs de G
288 //------------------------
289 // Ces coefficients sont egaux aux coefficients pairs du developpement de
290 // f.
291 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
292 // donnait exactement les coef. des cosinus, ce facteur serait 1.
293 
294  g[0] = cf0[0] ;
295  for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * cf0[i] ;
296  g[nm1] = cf0[nm1] ;
297 
298 // Transformation de Fourier inverse de G
299 //---------------------------------------
300 
301 // FFT inverse
302  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
303 
304 // Valeurs de f deduites de celles de G
305 //-------------------------------------
306 
307  for ( i = 1; i < nm1s2 ; i++ ) {
308 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
309  int isym = nm1 - i ;
310 // ... indice (dans le tableau ff0) du point x correspondant a psi
311  int ix = nm1 - i ;
312 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
313  int ixsym = nm1 - isym ;
314 
315  double fp = .5 * ( g[i] + g[isym] ) ;
316  double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
317 
318  ff0[ix] = fp + fm ;
319  ff0[ixsym] = fp - fm ;
320  }
321 
322 //... cas particuliers:
323  ff0[0] = g[0] - fmoins0 ;
324  ff0[nm1] = g[0] + fmoins0 ;
325  ff0[nm1s2] = g[nm1s2] ;
326 
327  } // fin de la boucle sur theta
328 
329 //--------------------------------------------------------------------
330 // partie sin(m phi) avec m pair : developpement en T_{2i}(x)
331 //--------------------------------------------------------------------
332 
333  j++ ;
334 
335  if ( (j != 1) && (j != n1f-1) ) {
336 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
337 // pas nuls
338 
339  for (k=0; k<n2c; k++) {
340 
341  int i0 = n2n3c * j + n3c * k ; // indice de depart
342  double* cf0 = cf + i0 ; // tableau des donnees a transformer
343 
344  i0 = n2n3f * j + n3f * k ; // indice de depart
345  double* ff0 = ff + i0 ; // tableau resultat
346 
347 /*
348  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
349  * reliee a x par x = cos(psi/2) et F(psi) = f(x(psi)).
350  */
351 
352 // Calcul des coefficients de Fourier de la fonction
353 // G(psi) = F+(psi) + F_(psi) sin(psi)
354 // en fonction des coefficients de Tchebyshev de f:
355 
356 // Coefficients impairs de G
357 //--------------------------
358 
359  double c1 = cf0[1] ;
360 
361  double som = 0;
362  ff0[1] = 0 ;
363  for ( i = 3; i < nr; i += 2 ) {
364  ff0[i] = cf0[i] - c1 ;
365  som += ff0[i] ;
366  }
367 
368 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
369  double fmoins0 = nm1s2 * c1 + som ;
370 
371 // Coef. impairs de G
372 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
373 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
374  g[1] = 0 ;
375  for ( i = 3; i < nr; i += 2 ) {
376  g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
377  }
378  g[nr] = 0 ;
379 
380 
381 // Coefficients pairs de G
382 //------------------------
383 // Ces coefficients sont egaux aux coefficients pairs du developpement de
384 // f.
385 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
386 // donnait exactement les coef. des cosinus, ce facteur serait 1.
387 
388  g[0] = cf0[0] ;
389  for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * cf0[i] ;
390  g[nm1] = cf0[nm1] ;
391 
392 // Transformation de Fourier inverse de G
393 //---------------------------------------
394 
395 // FFT inverse
396  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
397 
398 // Valeurs de f deduites de celles de G
399 //-------------------------------------
400 
401  for ( i = 1; i < nm1s2 ; i++ ) {
402 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
403  int isym = nm1 - i ;
404 // ... indice (dans le tableau ff0) du point x correspondant a psi
405  int ix = nm1 - i ;
406 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
407  int ixsym = nm1 - isym ;
408 
409  double fp = .5 * ( g[i] + g[isym] ) ;
410  double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
411 
412  ff0[ix] = fp + fm ;
413  ff0[ixsym] = fp - fm ;
414  }
415 
416 //... cas particuliers:
417  ff0[0] = g[0] - fmoins0 ;
418  ff0[nm1] = g[0] + fmoins0 ;
419  ff0[nm1s2] = g[nm1s2] ;
420 
421  } // fin de la boucle sur theta
422 
423  } // fin du cas ou le calcul etait necessaire (i.e. ou les
424  // coef en phi n'etaient pas nuls)
425 
426 // On passe au cas m pair suivant:
427  j+=3 ;
428 
429  } // fin de la boucle sur les cas m pair
430 
431  if (n1f<=3) { // cas m=0 seulement (symetrie axiale)
432  free (t1) ;
433  free (g) ;
434  return ;
435  }
436 
437 //=======================================================================
438 // Cas m impair
439 //=======================================================================
440 
441  j = 2 ;
442 
443  while (j<n1f-1) { //le dernier coef en phi n'est pas considere
444  // (car nul)
445 
446 //------------------------------------------------------------------------
447 // partie cos(m phi) avec m impair : developpement en T_{2i+1}(x)
448 //------------------------------------------------------------------------
449 
450  for (k=0; k<n2c; k++) {
451 
452  int i0 = n2n3c * j + n3c * k ; // indice de depart
453  double* cf0 = cf + i0 ; // tableau des donnees a transformer
454 
455  i0 = n2n3f * j + n3f * k ; // indice de depart
456  double* ff0 = ff + i0 ; // tableau resultat
457 
458 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
459 // h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le
460 // tableau t1 :
461  t1[0] = .5 * cf0[0] ;
462  for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ;
463  t1[nm1] = .5 * cf0[nr-2] ;
464 
465 /*
466  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
467  * reliee a x par x = cos(psi/2) et F(psi) = h(x(psi)).
468  */
469 
470 // Calcul des coefficients de Fourier de la fonction
471 // G(psi) = F+(psi) + F_(psi) sin(psi)
472 // en fonction des coefficients de Tchebyshev de f:
473 
474 // Coefficients impairs de G
475 //--------------------------
476 
477  double c1 = t1[1] ;
478 
479  double som = 0;
480  ff0[1] = 0 ;
481  for ( i = 3; i < nr; i += 2 ) {
482  ff0[i] = t1[i] - c1 ;
483  som += ff0[i] ;
484  }
485 
486 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
487  double fmoins0 = nm1s2 * c1 + som ;
488 
489 // Coef. impairs de G
490 // NB: le facteur 0.25 est du a la normalisation de fft991; si fft991
491 // donnait exactement les coef. des sinus, ce facteur serait -0.5.
492  g[1] = 0 ;
493  for ( i = 3; i < nr; i += 2 ) {
494  g[i] = 0.25 * ( ff0[i] - ff0[i-2] ) ;
495  }
496  g[nr] = 0 ;
497 
498 
499 // Coefficients pairs de G
500 //------------------------
501 // Ces coefficients sont egaux aux coefficients pairs du developpement de
502 // f.
503 // NB: le facteur 0.5 est du a la normalisation de fft991; si fft991
504 // donnait exactement les coef. des cosinus, ce facteur serait 1.
505 
506  g[0] = t1[0] ;
507  for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * t1[i] ;
508  g[nm1] = t1[nm1] ;
509 
510 // Transformation de Fourier inverse de G
511 //---------------------------------------
512 
513 // FFT inverse
514  F77_fft991( g, t1, trigo, facto, &inc, &jump, &nm1, &lot, &isign) ;
515 
516 // Valeurs de f deduites de celles de G
517 //-------------------------------------
518 
519  for ( i = 1; i < nm1s2 ; i++ ) {
520 // ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
521  int isym = nm1 - i ;
522 // ... indice (dans le tableau ff0) du point x correspondant a psi
523  int ix = nm1 - i ;
524 // ... indice (dans le tableau ff0) du point x correspondant a sym(psi)
525  int ixsym = nm1 - isym ;
526 
527  double fp = .5 * ( g[i] + g[isym] ) ;
528  double fm = .5 * ( g[i] - g[isym] ) / sinp[i] ;
529 
530  ff0[ix] = ( fp + fm ) / x[ix];
531  ff0[ixsym] = ( fp - fm ) / x[ixsym] ;
532  }
533 
534 //... cas particuliers:
535  ff0[0] = 0 ;
536  ff0[nm1] = g[0] + fmoins0 ;
537  ff0[nm1s2] = g[nm1s2] / x[nm1s2] ;
538 
539  } // fin de la boucle sur theta
540 
541 //------------------------------------------------------------------------
542 // partie sin(m phi) avec m impair : developpement en T_{2i+1}(x)
543 //------------------------------------------------------------------------
544 
545  j++ ;
546 
547  if ( j != n1f-1 ) {
548 // on effectue le calcul seulement dans les cas ou les coef en phi ne sont
549 // pas nuls
550 
551  for (k=0; k<n2c; k++) {
552 
553  int i0 = n2n3c * j + n3c * k ; // indice de depart
554  double* cf0 = cf + i0 ; // tableau des donnees a transformer
555 
556  i0 = n2n3f * j + n3f * k ; // indice de depart
557  double* ff0 = ff + i0 ; // tableau resultat
558 
559 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
560 // h(x) := x f(x) a partir des coefficients de f (resultat stoke dans le
561 // tableau t1 :
562  t1[0] = .5 * cf0[0] ;
563  for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[i] + cf0[i-1] ) ;
564  t1[nm1] = .5 * cf0[nr-2] ;
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) = h(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 = t1[1] ;
579 
580  double som = 0;
581  ff0[1] = 0 ;
582  for ( i = 3; i < nr; i += 2 ) {
583  ff0[i] = t1[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] = t1[0] ;
608  for (i=2; i<nm1; i += 2 ) g[i] = 0.5 * t1[i] ;
609  g[nm1] = t1[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 ) / x[ix];
632  ff0[ixsym] = ( fp - fm ) / x[ixsym] ;
633  }
634 
635 //... cas particuliers:
636  ff0[0] = 0 ;
637  ff0[nm1] = g[0] + fmoins0 ;
638  ff0[nm1s2] = g[nm1s2] / x[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 }
Lorene prototypes.
Definition: app_hor.h:67