LORENE
FFTW3/cftcossinci.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 * Transformation en cos((2*l+1)*theta) ou sin(2*l*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 cos( (2 l+1) theta ) .
69 * pour m impair:
70 * f(theta) = som_{l=0}^{nt-2} c_l sin( 2 l 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_{nt-1} = 0.
81 * Pour m impair, c_0 = 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: cftcossinci.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
90 * $Log: cftcossinci.C,v $
91 * Revision 1.4 2016/12/05 16:18:05 j_novak
92 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
93 *
94 * Revision 1.3 2014/10/13 08:53:19 j_novak
95 * Lorene classes and functions now belong to the namespace Lorene.
96 *
97 * Revision 1.2 2014/10/06 15:18:48 j_novak
98 * Modified #include directives to use c++ syntax.
99 *
100 * Revision 1.1 2004/12/21 17:06:02 j_novak
101 * Added all files for using fftw3.
102 *
103 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
104 * Suppressed the directive #include <malloc.h> for malloc is defined
105 * in <stdlib.h>
106 *
107 * Revision 1.3 2002/10/16 14:36:46 j_novak
108 * Reorganization of #include instructions of standard C++, in order to
109 * use experimental version 3 of gcc.
110 *
111 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon
112 * Modification of declaration of Fortran 77 prototypes for
113 * a better portability (in particular on IBM AIX systems):
114 * All Fortran subroutine names are now written F77_* and are
115 * defined in the new file C++/Include/proto_f77.h.
116 *
117 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
118 * LORENE
119 *
120 * Revision 2.0 1999/02/22 15:47:40 hyc
121 * *** empty log message ***
122 *
123 *
124 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cftcossinci.C,v 1.4 2016/12/05 16:18:05 j_novak Exp $
125 *
126 */
127
128
129// headers du C
130#include <cstdlib>
131#include <fftw3.h>
132
133//Lorene prototypes
134#include "tbl.h"
135
136// Prototypage des sous-routines utilisees:
137namespace Lorene {
138fftw_plan prepare_fft(int, Tbl*&) ;
139double* cheb_ini(const int) ;
140double* chebimp_ini(const int ) ;
141//*****************************************************************************
142
143void cftcossinci(const int* deg, const int* dimf, double* ff, const int* dimc,
144 double* cf)
145{
146
147int i, j, k ;
148
149// Dimensions des tableaux ff et cf :
150 int n1f = dimf[0] ;
151 int n2f = dimf[1] ;
152 int n3f = dimf[2] ;
153 int n1c = dimc[0] ;
154 int n2c = dimc[1] ;
155 int n3c = dimc[2] ;
156
157// Nombre de degres de liberte en theta :
158 int nt = deg[1] ;
159
160// Tests de dimension:
161 if (nt > n2f) {
162 cout << "cftcossinci: nt > n2f : nt = " << nt << " , n2f = "
163 << n2f << endl ;
164 abort () ;
165 }
166 if (nt > n2c) {
167 cout << "cftcossinci: nt > n2c : nt = " << nt << " , n2c = "
168 << n2c << endl ;
169 abort () ;
170 }
171 if (n1f > n1c) {
172 cout << "cftcossinci: n1f > n1c : n1f = " << n1f << " , n1c = "
173 << n1c << endl ;
174 abort () ;
175 }
176 if (n3f > n3c) {
177 cout << "cftcossinci: n3f > n3c : n3f = " << n3f << " , n3c = "
178 << n3c << endl ;
179 abort () ;
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 cos((2 l+1) 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// Multiplication de la fonction par x=cos(theta) (pour la rendre paire)
225// NB: dans les commentaires qui suivent, on note h(x) = x f(x).
226// (Les valeurs de h dans l'espace des configurations sont stokees dans le
227// tableau cf0).
228 for (i=0; i<nt-1; i++) cf0[n3c*i] = x[nm1-i] * ff0[n3f*i] ;
229 cf0[n3c*nm1] = 0 ;
230
231/*
232 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
233 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
234 */
235
236// Valeur en psi=0 de la partie antisymetrique de F, F_ :
237 double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
238
239// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
240//---------------------------------------------
241 for ( i = 1; i < nm1s2 ; i++ ) {
242// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
243 int isym = nm1 - i ;
244// ... indice (dans le tableau cf0) du point theta correspondant a psi
245 int ix = n3c * i ;
246// ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
247 int ixsym = n3c * isym ;
248// ... F+(psi)
249 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
250// ... F_(psi) sin(psi)
251 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
252 g.set(i) = fp + fms ;
253 g.set(isym) = fp - fms ;
254 }
255//... cas particuliers:
256 g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
257 g.set(nm1s2) = cf0[ n3c*nm1s2 ];
258
259// Developpement de G en series de Fourier par une FFT
260//----------------------------------------------------
261
262 fftw_execute(p) ;
263
264// Coefficients pairs du developmt. de Tchebyshev de h
265//----------------------------------------------------
266// Ces coefficients sont egaux aux coefficients en cosinus du developpement
267// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
268// de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
269// remplacer par un +1.) :
270
271 double fac = 2./double(nm1) ;
272 cf0[0] = g(0) / double(nm1) ;
273 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;
274 cf0[n3c*nm1] = g(nm1s2) ;
275
276// Coefficients impairs du developmt. de Tchebyshev de h
277//------------------------------------------------------
278// 1. Coef. c'_k (recurrence amorcee a partir de zero)
279// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
280// (si fftw donnait reellement les coef. en sinus, il faudrait le
281// remplacer par un -2.)
282 fac *= 2. ;
283 cf0[n3c] = 0 ;
284 double som = 0;
285 for ( i = 3; i < nt; i += 2 ) {
286 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
287 som += cf0[n3c*i] ;
288 }
289
290// 2. Calcul de c_1 :
291 double c1 = ( fmoins0 - som ) / nm1s2 ;
292
293// 3. Coef. c_k avec k impair:
294 cf0[n3c] = c1 ;
295 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
296
297
298// Coefficients de f en fonction de ceux de h
299//-------------------------------------------
300
301 cf0[0] = 2* cf0[0] ;
302 for (i=1; i<nm1; i++) {
303 cf0[n3c*i] = 2 * cf0[n3c*i] - cf0[n3c*(i-1)] ;
304 }
305 cf0[n3c*nm1] = 0 ;
306
307 } // fin de la boucle sur r
308
309
310//--------------------------------------------------------------------
311// partie sin(m phi) avec m pair : transformation en cos((2 l+1) theta)
312//--------------------------------------------------------------------
313
314 j++ ;
315
316 if ( (j != 1) && (j != n1f-1 ) ) {
317// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
318// pas nuls
319
320 for (k=0; k<n3f; k++) {
321
322 int i0 = n2n3f * j + k ; // indice de depart
323 double* ff0 = ff + i0 ; // tableau des donnees a transformer
324
325 i0 = n2n3c * j + k ; // indice de depart
326 double* cf0 = cf + i0 ; // tableau resultat
327
328// Multiplication de la fonction par x=cos(theta) (pour la rendre paire)
329// NB: dans les commentaires qui suivent, on note h(x) = x f(x).
330// (Les valeurs de h dans l'espace des configurations sont stokees dans le
331// tableau cf0).
332 for (i=0; i<nt-1; i++) cf0[n3c*i] = x[nm1-i] * ff0[n3f*i] ;
333 cf0[n3c*nm1] = 0 ;
334
335/*
336 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
337 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
338 */
339
340// Valeur en psi=0 de la partie antisymetrique de F, F_ :
341 double fmoins0 = 0.5 * ( cf0[0] - cf0[ n3c*nm1 ] );
342
343// Fonction G(psi) = F+(psi) + F_(psi) sin(psi)
344//---------------------------------------------
345 for ( i = 1; i < nm1s2 ; i++ ) {
346// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
347 int isym = nm1 - i ;
348// ... indice (dans le tableau cf0) du point theta correspondant a psi
349 int ix = n3c * i ;
350// ... indice (dans le tableau cf0) du point theta correspondant a sym(psi)
351 int ixsym = n3c * isym ;
352// ... F+(psi)
353 double fp = 0.5 * ( cf0[ix] + cf0[ixsym] ) ;
354// ... F_(psi) sin(psi)
355 double fms = 0.5 * ( cf0[ix] - cf0[ixsym] ) * sinp[i] ;
356 g.set(i) = fp + fms ;
357 g.set(isym) = fp - fms ;
358 }
359//... cas particuliers:
360 g.set(0) = 0.5 * ( cf0[0] + cf0[ n3c*nm1 ] );
361 g.set(nm1s2) = cf0[ n3c*nm1s2 ];
362
363// Developpement de G en series de Fourier par une FFT
364//----------------------------------------------------
365
366 fftw_execute(p) ;
367
368// Coefficients pairs du developmt. de Tchebyshev de h
369//----------------------------------------------------
370// Ces coefficients sont egaux aux coefficients en cosinus du developpement
371// de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
372// de fftw; si fftw donnait reellement les coef. en cosinus, il faudrait le
373// remplacer par un +1.) :
374
375 double fac = 2./double(nm1) ;
376 cf0[0] = g(0) / double(nm1) ;
377 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = fac * g(i/2) ;
378 cf0[n3c*nm1] = g(nm1s2) ;
379
380// Coefficients impairs du developmt. de Tchebyshev de h
381//------------------------------------------------------
382// 1. Coef. c'_k (recurrence amorcee a partir de zero)
383// Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
384// (si fftw donnait reellement les coef. en sinus, il faudrait le
385// remplacer par un -2.)
386 fac *= 2. ;
387 cf0[n3c] = 0 ;
388 double som = 0;
389 for ( i = 3; i < nt; i += 2 ) {
390 cf0[n3c*i] = cf0[n3c*(i-2)] + fac * g(nm1 - i/2) ;
391 som += cf0[n3c*i] ;
392 }
393
394// 2. Calcul de c_1 :
395 double c1 = ( fmoins0 - som ) / nm1s2 ;
396
397// 3. Coef. c_k avec k impair:
398 cf0[n3c] = c1 ;
399 for ( i = 3; i < nt; i += 2 ) cf0[n3c*i] += c1 ;
400
401
402// Coefficients de f en fonction de ceux de h
403//-------------------------------------------
404
405 cf0[0] = 2* cf0[0] ;
406 for (i=1; i<nm1; i++) {
407 cf0[n3c*i] = 2 * cf0[n3c*i] - cf0[n3c*(i-1)] ;
408 }
409 cf0[n3c*nm1] = 0 ;
410
411 } // fin de la boucle sur r
412
413 } // fin du cas ou le calcul etait necessaire (i.e. ou les
414 // coef en phi n'etaient pas nuls)
415
416// On passe au cas m pair suivant:
417 j+=3 ;
418
419 } // fin de la boucle sur les cas m pair
420
421 if (n1f<=3) // cas m=0 seulement (symetrie axiale)
422 return ;
423
424//=======================================================================
425// Cas m impair
426//=======================================================================
427
428 j = 2 ;
429
430 while (j<n1f-1) { //le dernier coef en phi n'est pas considere
431 // (car nul)
432
433//--------------------------------------------------------------------
434// partie cos(m phi) avec m impair : transformation en sin((2 l) theta)
435//--------------------------------------------------------------------
436
437 for (k=0; k<n3f; k++) {
438
439 int i0 = n2n3f * j + k ; // indice de depart
440 double* ff0 = ff + i0 ; // tableau des donnees a transformer
441
442 i0 = n2n3c * j + k ; // indice de depart
443 double* cf0 = cf + i0 ; // tableau resultat
444
445/*
446 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
447 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
448 */
449
450// Fonction G(psi) = F+(psi) sin(psi) + F_(psi)
451//---------------------------------------------
452 for ( i = 1; i < nm1s2 ; i++ ) {
453// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
454 int isym = nm1 - i ;
455// ... indice (dans le tableau ff0) du point theta correspondant a psi
456 int ix = n3f * i ;
457// ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
458 int ixsym = n3f * isym ;
459// ... F+(psi) sin(psi)
460 double fps = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
461// ... F_(psi)
462 double fm = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
463 g.set(i) = fps + fm ;
464 g.set(isym) = fps - fm ;
465 }
466//... cas particuliers:
467 g.set(0) = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ) ;
468 g.set(nm1s2) = ff0[ n3f*nm1s2 ] ;
469
470// Developpement de G en series de Fourier par une FFT
471//----------------------------------------------------
472 fftw_execute(p) ;
473
474// Coefficients pairs du developmt. sin(2l theta) de f
475//----------------------------------------------------
476// Ces coefficients sont egaux aux coefficients en sinus du developpement
477// de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
478// de fftw: si fftw donnait reellement les coefficients en sinus,
479// il faudrait le remplacer par un +1) :
480 double fac = 2./double(nm1) ;
481 cf0[0] = 0. ;
482 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = -fac * g(nm1 - i/2);
483 cf0[n3c*nm1] = 0. ;
484
485// Coefficients impairs du developmt. en sin(2l theta) de f
486//---------------------------------------------------------
487// Ces coefficients sont obtenus par une recurrence a partir des
488// coefficients en cosinus du developpement de G en series de Fourier
489// (le facteur 4/nm1 vient de la normalisation
490// de fftw: si fftw donnait reellement les coefficients en cosinus,
491// il faudrait le remplacer par un +2.)
492
493 cf0[n3c] = fac * g(0) ;
494 fac *= 2. ;
495 for ( i = 3; i < nt; i += 2 ) {
496 cf0[ n3c*i ] = cf0[ n3c*(i-2) ] + fac * g(i/2) ;
497 }
498
499 } // fin de la boucle sur r
500
501//------------------------------------------------------------------------
502// partie sin(m phi) avec m impair : transformation en sin((2 l) theta)
503//------------------------------------------------------------------------
504
505 j++ ;
506
507 if ( j != n1f-1 ) {
508// on effectue le calcul seulement dans les cas ou les coef en phi ne sont
509// pas nuls
510
511 for (k=0; k<n3f; k++) {
512
513 int i0 = n2n3f * j + k ; // indice de depart
514 double* ff0 = ff + i0 ; // tableau des donnees a transformer
515
516 i0 = n2n3c * j + k ; // indice de depart
517 double* cf0 = cf + i0 ; // tableau resultat
518
519/*
520 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
521 * reliee a theta par psi = 2 theta et F(psi) = f(theta(psi)).
522 */
523
524// Fonction G(psi) = F+(psi) sin(psi) + F_(psi)
525//---------------------------------------------
526 for ( i = 1; i < nm1s2 ; i++ ) {
527// ... indice (dans le tableau g) du pt symetrique de psi par rapport a pi/2:
528 int isym = nm1 - i ;
529// ... indice (dans le tableau ff0) du point theta correspondant a psi
530 int ix = n3f * i ;
531// ... indice (dans le tableau ff0) du point theta correspondant a sym(psi)
532 int ixsym = n3f * isym ;
533// ... F+(psi) sin(psi)
534 double fps = 0.5 * ( ff0[ix] + ff0[ixsym] ) * sinp[i] ;
535// ... F_(psi)
536 double fm = 0.5 * ( ff0[ix] - ff0[ixsym] ) ;
537 g.set(i) = fps + fm ;
538 g.set(isym) = fps - fm ;
539 }
540//... cas particuliers:
541 g.set(0) = 0.5 * ( ff0[0] - ff0[ n3f*nm1 ] ) ;
542 g.set(nm1s2) = ff0[ n3f*nm1s2 ] ;
543
544// Developpement de G en series de Fourier par une FFT
545//----------------------------------------------------
546
547 fftw_execute(p) ;
548
549// Coefficients pairs du developmt. sin(2l theta) de f
550//----------------------------------------------------
551// Ces coefficients sont egaux aux coefficients en sinus du developpement
552// de G en series de Fourier (le facteur -2/nm1 vient de la normalisation
553// de fftw: si fftw donnait reellement les coefficients en sinus,
554// il faudrait le remplacer par un +1) :
555 double fac = 2./double(nm1) ;
556 cf0[0] = 0. ;
557 for (i=2; i<nm1; i += 2 ) cf0[n3c*i] = -fac * g(nm1 - i/2);
558 cf0[n3c*nm1] = 0. ;
559
560// Coefficients impairs du developmt. en sin(2l theta) de f
561//---------------------------------------------------------
562// Ces coefficients sont obtenus par une recurrence a partir des
563// coefficients en cosinus du developpement de G en series de Fourier
564// (le facteur 4/nm1 vient de la normalisation
565// de fftw: si fftw donnait reellement les coefficients en cosinus,
566// il faudrait le remplacer par un +2.)
567
568 cf0[n3c] = fac * g(0) ;
569 fac *= 2. ;
570 for ( i = 3; i < nt; i += 2 ) {
571 cf0[ n3c*i ] = cf0[ n3c*(i-2) ] + fac * g(i/2) ;
572 }
573
574 } // fin de la boucle sur r
575
576 } // fin du cas ou le calcul etait necessaire (i.e. ou les
577 // coef en phi n'etaient pas nuls)
578
579
580// On passe au cas m impair suivant:
581 j+=3 ;
582
583 } // fin de la boucle sur les cas m impair
584
585
586}
587}
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67
Coord x
x coordinate centered on the grid
Definition map.h:738
Coord sinp
Definition map.h:735