LORENE
scalar_manip.C
1 /*
2  * Copyright (c) 2003 Eric Gourgoulhon & Jerome Novak
3  *
4  * Copyright (c) 2000-2001 Philippe Grandclement (Cmp version)
5  *
6  * This file is part of LORENE.
7  *
8  * LORENE is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * LORENE is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with LORENE; if not, write to the Free Software
20  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
21  *
22  */
23 
24 
25 
26 
27 /*
28  * $Id: scalar_manip.C,v 1.21 2018/11/16 14:34:37 j_novak Exp $
29  * $Log: scalar_manip.C,v $
30  * Revision 1.21 2018/11/16 14:34:37 j_novak
31  * Changed minor points to avoid some compilation warnings.
32  *
33  * Revision 1.20 2016/12/05 16:18:18 j_novak
34  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35  *
36  * Revision 1.19 2014/10/13 08:53:46 j_novak
37  * Lorene classes and functions now belong to the namespace Lorene.
38  *
39  * Revision 1.18 2014/10/06 15:16:15 j_novak
40  * Modified #include directives to use c++ syntax.
41  *
42  * Revision 1.17 2008/10/03 09:03:52 j_novak
43  * Correction of yet another mistake (the array values in physical space was not
44  * destroyed).
45  *
46  * Revision 1.16 2008/09/30 08:35:18 j_novak
47  * Correction of forgotten call to coef()
48  *
49  * Revision 1.15 2008/09/29 13:23:51 j_novak
50  * Implementation of the angular mapping associated with an affine
51  * mapping. Things must be improved to take into account the domain index.
52  *
53  * Revision 1.14 2008/09/22 19:08:01 j_novak
54  * New methods to deal with boundary conditions
55  *
56  * Revision 1.13 2007/06/21 20:00:00 k_taniguchi
57  * Addition of the method filtre_r (int n, int nz)
58  *
59  * Revision 1.12 2006/06/28 07:46:41 j_novak
60  * Better treatment in the case of a domain set to zero.
61  *
62  * Revision 1.11 2005/09/07 13:39:10 j_novak
63  * *** empty log message ***
64  *
65  * Revision 1.10 2005/09/07 13:10:48 j_novak
66  * Added a filter setting to zero all mulitpoles in a given range.
67  *
68  * Revision 1.9 2004/11/23 12:47:44 f_limousin
69  * Add function filtre_tp(int nn, int nz1, int nz2).
70  *
71  * Revision 1.8 2004/05/07 11:24:54 f_limousin
72  * Implement new method filtre_r (int* nn)
73  *
74  * Revision 1.7 2004/02/27 09:47:26 f_limousin
75  * New methods filtre_phi(int) and filtre_theta(int).
76  *
77  * Revision 1.6 2004/01/23 13:26:28 e_gourgoulhon
78  * Added methods set_inner_boundary and set_outer_boundary.
79  * Methods set_val_inf and set_val_hor, which are particular cases of
80  * the above, have been suppressed.
81  *
82  * Revision 1.5 2003/11/13 13:43:55 p_grandclement
83  * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
84  *
85  * Revision 1.4 2003/10/11 14:47:17 e_gourgoulhon
86  * Lines 112 and 145 : replaced "0" by "double(0)" in the comparison
87  * statement.
88  *
89  * Revision 1.3 2003/10/10 15:57:29 j_novak
90  * Added the state one (ETATUN) to the class Scalar
91  *
92  * Revision 1.2 2003/10/08 14:24:09 j_novak
93  * replaced mult_r_zec with mult_r_ced
94  *
95  * Revision 1.1 2003/09/25 09:33:36 j_novak
96  * Added methods for integral calculation and various manipulations
97  *
98  *
99  * $Header: /cvsroot/Lorene/C++/Source/Tensor/Scalar/scalar_manip.C,v 1.21 2018/11/16 14:34:37 j_novak Exp $
100  *
101  */
102 
103 //standard
104 #include <cstdlib>
105 #include <cmath>
106 
107 // Lorene
108 #include "tensor.h"
109 #include "proto.h"
110 #include "utilitaires.h"
111 
112 /*
113  * Annule tous les l entre l_min et l_max (compris)
114  */
115 
116 namespace Lorene {
117 void Scalar::annule_l (int l_min, int l_max, bool ylm_output) {
118 
119  assert (etat != ETATNONDEF) ;
120  assert (l_min <= l_max) ;
121  assert (l_min >= 0) ;
122  if (etat == ETATZERO )
123  return ;
124 
125  if (etat == ETATUN) {
126  if (l_min == 0) set_etat_zero() ;
127  else return ;
128  }
129 
130  va.ylm() ;
131  Mtbl_cf& m_coef = *va.c_cf ;
132  const Base_val& base = va.base ;
133  int l_q, m_q, base_r ;
134  for (int lz=0; lz<mp->get_mg()->get_nzone(); lz++)
135  if (m_coef(lz).get_etat() != ETATZERO)
136  for (int k=0; k<mp->get_mg()->get_np(lz)+1; k++)
137  for (int j=0; j<mp->get_mg()->get_nt(lz); j++)
138  for (int i=0; i<mp->get_mg()->get_nr(lz); i++) {
139  base.give_quant_numbers(lz, k, j, m_q, l_q, base_r) ;
140  if ((l_min <= l_q) && (l_q<= l_max))
141  m_coef.set(lz, k, j, i) = 0 ;
142  }
143  if (va.c != 0x0) {
144  delete va.c ;
145  va.c = 0x0 ;
146  }
147  if (!ylm_output)
148  va.ylm_i() ;
149 
150  return ;
151 }
152 
153 /*
154  * Annule les n derniers coefficients en r dans la derniere zone
155  */
156 
157 void Scalar::filtre (int n) {
158 
159  assert (etat != ETATNONDEF) ;
160  if ( (etat == ETATZERO) || (etat == ETATUN) )
161  return ;
162 
163  int nz = mp->get_mg()->get_nzone() ;
164  int np = mp->get_mg()->get_np(nz-1) ;
165  int nt = mp->get_mg()->get_nt(nz-1) ;
166  int nr = mp->get_mg()->get_nr(nz-1) ;
167 
168  del_deriv() ;
169 
170  va.coef() ;
171  va.set_etat_cf_qcq() ;
172 
173 
174  for (int k=0 ; k<np+1 ; k++)
175  if (k!=1)
176  for (int j=0 ; j<nt ; j++)
177  for (int i=nr-1 ; i>nr-1-n ; i--)
178  va.c_cf->set(nz-1, k, j, i) = 0 ;
179 }
180 
181 
182 /*
183  * Annule les n derniers coefficients en r dans toutes les zones
184  */
185 
186 void Scalar::filtre_r (int* nn) {
187  assert (etat != ETATNONDEF) ;
188  if ( (etat == ETATZERO) || (etat == ETATUN) )
189  return ;
190 
191  del_deriv() ;
192 
193  va.coef() ;
194  va.set_etat_cf_qcq() ;
195  int nz = mp->get_mg()->get_nzone() ;
196  int* nr = new int[nz];
197  int* nt = new int[nz];
198  int* np = new int[nz];
199  for (int l=0; l<=nz-1; l++) {
200  nr[l] = mp->get_mg()->get_nr(l) ;
201  nt[l] = mp->get_mg()->get_nt(l) ;
202  np[l] = mp->get_mg()->get_np(l) ;
203  }
204 
205  for (int l=0; l<=nz-1; l++)
206  for (int k=0 ; k<np[l]+1 ; k++)
207  if (k!=1)
208  for (int j=0 ; j<nt[l] ; j++)
209  for (int i=nr[l]-1; i>nr[l]-1-nn[l] ; i--)
210  va.c_cf->set(l, k, j, i) = 0 ;
211 
212  if (va.c != 0x0) {
213  delete va.c ;
214  va.c = 0x0 ;
215  }
216 
217 }
218 
219 
220 /*
221  * Annule les n derniers coefficients en r dans zone nz
222  */
223 
224 void Scalar::filtre_r (int n, int nz) {
225  assert (etat != ETATNONDEF) ;
226  if ( (etat == ETATZERO) || (etat == ETATUN) )
227  return ;
228 
229  del_deriv() ;
230 
231  va.coef() ;
232  va.set_etat_cf_qcq() ;
233  int nr = mp->get_mg()->get_nr(nz) ;
234  int nt = mp->get_mg()->get_nt(nz) ;
235  int np = mp->get_mg()->get_np(nz) ;
236 
237  for (int k=0 ; k<np+1 ; k++)
238  if (k!=1)
239  for (int j=0 ; j<nt ; j++)
240  for (int i=nr-1; i>nr-1-n ; i--)
241  va.c_cf->set(nz, k, j, i) = 0 ;
242 
243  if (va.c != 0x0) {
244  delete va.c ;
245  va.c = 0x0 ;
246  }
247 
248 }
249 
250 
251 /*
252  * Annule les n derniers coefficients en phi dans zone nz
253  */
254 
255 void Scalar::filtre_phi (int n, int nz) {
256  assert (etat != ETATNONDEF) ;
257  if ( (etat == ETATZERO) || (etat == ETATUN) )
258  return ;
259 
260  del_deriv() ;
261 
262  va.coef() ;
263  va.set_etat_cf_qcq() ;
264  int np = mp->get_mg()->get_np(nz) ;
265  int nt = mp->get_mg()->get_nt(nz) ;
266  int nr = mp->get_mg()->get_nr(nz) ;
267 
268  for (int k=np+1-n ; k<np+1 ; k++)
269  for (int j=0 ; j<nt ; j++)
270  for (int i=0 ; i<nr ; i++)
271  va.c_cf->set(nz, k, j, i) = 0 ;
272 
273 }
274 
275 
276 void Scalar::filtre_tp(int nn, int nz1, int nz2) {
277 
278  va.filtre_tp(nn, nz1, nz2) ;
279 
280 }
281 
282 
283  /* Sets the value of the {\tt Scalar} at the inner boundary of a given
284  * domain.
285  * @param l [input] domain index
286  * @param x [input] (constant) value at the inner boundary of domain no. {\tt l}
287  */
288 
289 void Scalar::set_inner_boundary(int l0, double x0) {
290 
291  assert (etat != ETATNONDEF) ;
292  if (etat == ETATZERO) {
293  if (x0 == double(0)) return ;
294  else annule_hard() ;
295  }
296 
297  if (etat == ETATUN) {
298  if (x0 == double(1)) return ;
299  else etat = ETATQCQ ;
300  }
301 
302  del_deriv() ;
303 
304  int nt = mp->get_mg()->get_nt(l0) ;
305  int np = mp->get_mg()->get_np(l0) ;
306 
307  va.coef_i() ;
308  va.set_etat_c_qcq() ;
309 
310  for (int k=0 ; k<np ; k++)
311  for (int j=0 ; j<nt ; j++)
312  va.set(l0, k, j, 0) = x0 ;
313 }
314 
315  /* Sets the value of the {\tt Scalar} at the outer boundary of a given
316  * domain.
317  * @param l [input] domain index
318  * @param x [input] (constant) value at the outer boundary of domain no. {\tt l}
319  */
320 
321 void Scalar::set_outer_boundary(int l0, double x0) {
322 
323  assert (etat != ETATNONDEF) ;
324  if (etat == ETATZERO) {
325  if (x0 == double(0)) return ;
326  else annule_hard() ;
327  }
328 
329  if (etat == ETATUN) {
330  if (x0 == double(1)) return ;
331  else etat = ETATQCQ ;
332  }
333 
334  del_deriv() ;
335 
336  int nrm1 = mp->get_mg()->get_nr(l0) - 1 ;
337  int nt = mp->get_mg()->get_nt(l0) ;
338  int np = mp->get_mg()->get_np(l0) ;
339 
340  va.coef_i() ;
341  va.set_etat_c_qcq() ;
342 
343  for (int k=0 ; k<np ; k++)
344  for (int j=0 ; j<nt ; j++)
345  va.set(l0, k, j, nrm1) = x0 ;
346 }
347 
348 /*
349  * Permet de fixer la decroissance du cmp a l infini en viurant les
350  * termes en 1/r^n
351  */
352 void Scalar::fixe_decroissance (int puis) {
353 
354  if (puis<dzpuis)
355  return ;
356  else {
357 
358  int nbre = puis-dzpuis ;
359 
360  // le confort avant tout ! (c'est bien le confort ...)
361  int nz = mp->get_mg()->get_nzone() ;
362  int np = mp->get_mg()->get_np(nz-1) ;
363  int nt = mp->get_mg()->get_nt(nz-1) ;
364  int nr = mp->get_mg()->get_nr(nz-1) ;
365 
366  const Map_af* map = dynamic_cast<const Map_af*>(mp) ;
367  if (map == 0x0) {
368  cout << "Le mapping doit etre affine" << endl ;
369  abort() ;
370  }
371 
372  double alpha = map->get_alpha()[nz-1] ;
373 
374  Scalar courant (*this) ;
375 
376  va.coef() ;
377  va.set_etat_cf_qcq() ;
378 
379  for (int conte=0 ; conte<nbre ; conte++) {
380 
381  int base_r = courant.va.base.get_base_r(nz-1) ;
382 
383  courant.va.coef() ;
384 
385  // On calcul les coefficients de 1/r^conte
386  double* coloc = new double [nr] ;
387  int * deg = new int[3] ;
388  deg[0] = 1 ;
389  deg[1] = 1 ;
390  deg[2] = nr ;
391 
392  for (int i=0 ; i<nr ; i++)
393  coloc[i] =pow(alpha, double(conte))*
394  pow(-1-cos(M_PI*i/(nr-1)), double(conte)) ;
395 
396  cfrcheb(deg, deg, coloc, deg, coloc) ;
397 
398  for (int k=0 ; k<np+1 ; k++)
399  if (k != 1)
400  for (int j=0 ; j<nt ; j++) {
401 
402  // On doit determiner le coefficient du truc courant :
403  double* coef = new double [nr] ;
404  double* auxi = new double[1] ;
405  for (int i=0 ; i<nr ; i++)
406  coef[i] = (*courant.va.c_cf)(nz-1, k, j, i) ;
407  switch (base_r) {
408  case R_CHEBU :
409  som_r_chebu (coef, nr, 1, 1, 1, auxi) ;
410  break ;
411  default :
412  som_r_pas_prevu (coef, nr, 1, 1, 1, auxi) ;
413  break ;
414  }
415 
416  // On modifie le cmp courant :
417  courant.va.coef() ;
418  courant.va.set_etat_cf_qcq() ;
419  courant.va.c_cf->set(nz-1, k, j, 0) -= *auxi ;
420 
421  for (int i=0 ; i<nr ; i++)
422  this->va.c_cf->set(nz-1, k, j, i) -= *auxi * coloc[i] ;
423 
424 
425  delete [] coef ;
426  delete [] auxi ;
427  }
428  delete [] coloc ;
429  delete [] deg ;
430 
431  courant.mult_r_ced() ;
432  }
433  }
434 }
435 
436 Tbl Scalar::tbl_out_bound(int l_zone, bool output_ylm) {
437  va.coef() ;
438  if (output_ylm) va.ylm() ;
439 
440  int np = mp->get_mg()->get_np(l_zone) ;
441  int nt = mp->get_mg()->get_nt(l_zone) ;
442  Tbl resu(np+2, nt) ;
443  if (etat == ETATZERO) resu.set_etat_zero() ;
444  else {
445  assert(etat == ETATQCQ) ;
446  resu.set_etat_qcq() ;
447  for (int k=0; k<np+2; k++)
448  for (int j=0; j<nt; j++)
449  resu.set(k, j) = va.c_cf->val_out_bound_jk(l_zone, j, k) ;
450  }
451  return resu ;
452 }
453 
454 Tbl Scalar::tbl_in_bound(int l_zone, bool output_ylm) {
455  assert(mp->get_mg()->get_type_r(l_zone) != RARE) ;
456  va.coef() ;
457  if (output_ylm) va.ylm() ;
458 
459  int np = mp->get_mg()->get_np(l_zone) ;
460  int nt = mp->get_mg()->get_nt(l_zone) ;
461  Tbl resu(np+2, nt) ;
462  if (etat == ETATZERO) resu.set_etat_zero() ;
463  else {
464  assert(etat == ETATQCQ) ;
465  resu.set_etat_qcq() ;
466  for (int k=0; k<np+2; k++)
467  for (int j=0; j<nt; j++)
468  resu.set(k, j) = va.c_cf->val_in_bound_jk(l_zone, j, k) ;
469  }
470  return resu ;
471 }
472 
473 Scalar Scalar::scalar_out_bound(int l_zone, bool output_ylm) {
474  va.coef() ;
475  if (output_ylm) va.ylm() ;
476 
477  Scalar resu(mp->mp_angu(l_zone)) ;
478  resu.std_spectral_base() ;
479  Base_val base = resu.get_spectral_base() ;
480  base.set_base_t(va.base.get_base_t(l_zone)) ;
481  resu.set_spectral_base(base) ;
482  if (etat == ETATZERO) resu.set_etat_zero() ;
483  else {
484  assert(etat == ETATQCQ) ;
485  resu.annule_hard() ;
486  int np = mp->get_mg()->get_np(l_zone) ;
487  int nt = mp->get_mg()->get_nt(l_zone) ;
488  for (int k=0; k<np+2; k++)
489  for (int j=0; j<nt; j++)
490  resu.set_spectral_va().c_cf->set(0, k, j, 0)
491  = va.c_cf->val_out_bound_jk(l_zone, j, k) ;
492  delete resu.set_spectral_va().c ;
493  resu.set_spectral_va().c = 0x0 ;
494  }
495  return resu ;
496 }
497 }
void annule_l(int l_min, int l_max, bool ylm_output=false)
Sets all the multipolar components between l_min and l_max to zero.
Definition: scalar_manip.C:117
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
void filtre_r(int *nn)
Sets the n lasts coefficients in r to 0 in all domains.
Definition: scalar_manip.C:186
const double * get_alpha() const
Returns the pointer on the array alpha.
Definition: map_af.C:604
virtual const Map_af & mp_angu(int) const =0
Returns the "angular" mapping for the outside of domain l_zone.
void ylm_i()
Inverse of ylm()
Definition: valeur_ylm_i.C:134
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
void set_etat_cf_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_c...
Definition: valeur.C:715
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:151
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
void give_quant_numbers(int, int, int, int &, int &, int &) const
Computes the various quantum numbers and 1d radial base.
Lorene prototypes.
Definition: app_hor.h:67
Tbl & set(int l)
Read/write of the Tbl containing the coefficients in a given domain.
Definition: mtbl_cf.h:304
int get_base_t(int l) const
Returns the expansion basis for functions in the domain of index l (e.g.
Definition: base_val.h:414
void ylm()
Computes the coefficients of *this.
Definition: valeur_ylm.C:141
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:777
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
void coef_i() const
Computes the physical value of *this.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
Scalar scalar_out_bound(int n, bool leave_ylm=false)
Returns the Scalar containing the values of angular coefficients at the outer boundary.
Definition: scalar_manip.C:473
int dzpuis
Power of r by which the quantity represented by this must be divided in the compactified external do...
Definition: scalar.h:409
void set_base_t(int base_t)
Sets the expansion basis for functions in all domains.
Definition: base_val.C:198
Tbl tbl_in_bound(int n, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the inner boundary.
Definition: scalar_manip.C:454
int get_etat() const
Returns the logical state.
Definition: mtbl_cf.h:466
double val_out_bound_jk(int l, int j, int k) const
Computes the angular coefficient of index j,k of the field represented by *this at by means of the s...
void annule_hard()
Sets the Scalar to zero in a hard way.
Definition: scalar.C:386
int get_base_r(int l) const
Returns the expansion basis for r ( ) functions in the domain of index l (e.g.
Definition: base_val.h:403
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
friend Scalar pow(const Scalar &, int)
Power .
Definition: scalar_math.C:454
double val_in_bound_jk(int l, int j, int k) const
Computes the angular coefficient of index j,k of the field represented by *this at by means of the s...
friend Scalar cos(const Scalar &)
Cosine.
Definition: scalar_math.C:107
virtual void del_deriv() const
Logical destructor of the derivatives.
Definition: scalar.C:293
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
Definition: scalar_manip.C:321
void filtre_tp(int nn, int nz1, int nz2)
Sets the n lasts coefficients in to 0 in the domain nz1 to nz2 when expressed in spherical harmonics...
Definition: valeur.C:927
void mult_r_ced()
Multiplication by r in the compactified external domain (CED), the dzpuis flag is not changed...
Mtbl * c
Values of the function at the points of the multi-grid.
Definition: valeur.h:309
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
Valeur va
The numerical value of the Scalar.
Definition: scalar.h:411
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition: scalar_manip.C:157
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: tbl.C:350
void set_inner_boundary(int l, double x)
Sets the value of the Scalar at the inner boundary of a given domain.
Definition: scalar_manip.C:289
void filtre_tp(int nn, int nz1, int nz2)
Sets the n lasts coefficients in to 0 in the domain nz1 to nz2 when expressed in spherical harmonics...
Definition: scalar_manip.C:276
void fixe_decroissance(int puis)
Substracts all the components behaving like in the external domain, with n strictly lower than puis ...
Definition: scalar_manip.C:352
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
Bases of the spectral expansions.
Definition: base_val.h:325
void filtre_phi(int n, int zone)
Sets the n lasts coefficients in to 0 in the domain zone .
Definition: scalar_manip.C:255
Affine radial mapping.
Definition: map.h:2042
int etat
The logical state ETATNONDEF (undefined), ETATZERO (null), ETATUN (one), or ETATQCQ (ordinary)...
Definition: scalar.h:402
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c )...
Definition: valeur.C:704
Tbl tbl_out_bound(int l_dom, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the outer boundary.
Definition: scalar_manip.C:436
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:474
#define R_CHEBU
base de Chebychev ordinaire (fin), dev. en 1/r
Definition: type_parite.h:180
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:491
const Map *const mp
Mapping on which the numerical values at the grid points are defined.
Definition: tensor.h:301
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition: valeur.h:373