LORENE
hot_star.C
1 /*
2  * Methods of class Hot_star
3  *
4  * (see file star.h for documentation)
5  */
6 
7 /*
8  * Copyright (c) 2004 Francois Limousin
9  *
10  * Copyright (c) 2000-2001 Eric Gourgoulhon (for preceding class Etoile)
11  * Copyright (c) 2000-2001 Keisuke Taniguchi (for preceding class Etoile)
12  *
13  * This file is part of LORENE.
14  *
15  * LORENE is free software; you can redistribute it and/or modify
16  * it under the terms of the GNU General Public License as published by
17  * the Free Software Foundation; either version 2 of the License, or
18  * (at your option) any later version.
19  *
20  * LORENE is distributed in the hope that it will be useful,
21  * but WITHOUT ANY WARRANTY; without even the implied warranty of
22  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23  * GNU General Public License for more details.
24  *
25  * You should have received a copy of the GNU General Public License
26  * along with LORENE; if not, write to the Free Software
27  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28  *
29  */
30 
31 
32 
33 
34 /*
35  * $Id: hot_star.C,v 1.1 2026/03/23 17:33:04 s_jaraba Exp $
36  * $Log: hot_star.C,v $
37  * Revision 1.1 2026/03/23 17:33:04 s_jaraba
38  * Initial commit for classes Hot_star, Hot_star_rot_CFC and Hot_star_rot_diff_CFC
39  *
40  * Revision 1.22 2016/12/05 16:18:14 j_novak
41  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
42  *
43  * Revision 1.21 2014/10/13 08:53:37 j_novak
44  * Lorene classes and functions now belong to the namespace Lorene.
45  *
46  * Revision 1.20 2013/04/20 20:56:15 m_bejger
47  * Fix for three domains in star in Star::equation_of_state from Etoile/etoile.C
48  *
49  * Revision 1.19 2010/02/02 12:45:16 e_gourgoulhon
50  * Improved the display (operator>>)
51  *
52  * Revision 1.18 2010/01/26 16:49:03 e_gourgoulhon
53  * Commented the test on the relativistic character of the EOS: the
54  * relativity parameter is not defined (yet !) in the base class Star.
55  *
56  * Revision 1.17 2007/11/06 16:22:03 j_novak
57  * The data member stress_euler is now a Sym_tensor instead of a Tensor.
58  *
59  * Revision 1.16 2007/06/21 19:53:47 k_taniguchi
60  * Addition of p_ray_eq_3pis2
61  *
62  * Revision 1.15 2005/09/14 12:30:52 f_limousin
63  * Saving of fields lnq and logn in class Star.
64  *
65  * Revision 1.14 2005/09/13 19:38:31 f_limousin
66  * Reintroduction of the resolution of the equations in cartesian coordinates.
67  *
68  * Revision 1.13 2005/02/17 17:29:04 f_limousin
69  * Change the name of some quantities to be consistent with other classes
70  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
71  *
72  * Revision 1.12 2005/01/05 17:43:03 f_limousin
73  * u_euler is now constructed in the spherical triad to be consistent
74  * with all the others vectors ans tensors.
75  *
76  * Revision 1.11 2004/11/11 16:29:49 j_novak
77  * set_der_0x0 is no longer virtual (to be coherent with Tensor/Scalar classes).
78  *
79  * Revision 1.10 2004/06/22 12:48:08 f_limousin
80  * Change qq, qq_auto and qq_comp to beta, beta_auto and beta_comp.
81  *
82  * Revision 1.9 2004/06/07 16:21:35 f_limousin
83  * Add outputs
84  *
85  * Revision 1.8 2004/04/08 16:32:10 f_limousin
86  * The new variable is ln(Q) instead of Q=psi^2*N. It improves the
87  * convergence of the code.
88  *
89  * Revision 1.7 2004/03/25 10:29:26 j_novak
90  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
91  *
92  * Revision 1.6 2004/03/08 11:48:00 f_limousin
93  * Error in del_deriv() and set_der_0x0() : p_mass_b and p_mass_g were
94  * missing. And so they were never recomputed.
95  *
96  * Revision 1.5 2004/02/27 09:36:46 f_limousin
97  * u_euler is now constructed on a cartesian basis instead
98  * of a spherical one.
99  *
100  * Revision 1.4 2004/02/21 17:05:13 e_gourgoulhon
101  * Method Scalar::point renamed Scalar::val_grid_point.
102  * Method Scalar::set_point renamed Scalar::set_grid_point.
103  *
104  * Revision 1.3 2004/01/20 15:16:58 f_limousin
105  * First version
106  *
107  *
108  * $Header: /cvsroot/Lorene/C++/Source/Star/hot_star.C,v 1.1 2026/03/23 17:33:04 s_jaraba Exp $
109  *
110  */
111 
112 // Headers C
113 #include "math.h"
114 
115 // Headers Lorene
116 #include "hot_star.h"
117 #include "hoteos.h"
118 #include "utilitaires.h"
119 #include "param.h"
120 #include "unites.h"
121 
122 
123  //--------------//
124  // Constructors //
125  //--------------//
126 
127 // Standard constructor
128 // --------------------
129 namespace Lorene {
130 Hot_star::Hot_star(Map& mpi, int nzet_i, const Hot_eos& heos_i)
131  : mp(mpi),
132  nzet(nzet_i),
133  heos(heos_i),
134  ent(mpi),
135  nbar(mpi),
136  ener(mpi),
137  press(mpi),
138  temp(mpi),
139  entropy(mpi),
140  ener_euler(mpi),
141  s_euler(mpi),
142  gam_euler(mpi),
143  u_euler(mpi, CON, mp.get_bvect_spher()),
144  stress_euler(mpi, CON, mp.get_bvect_spher()),
145  logn(mpi),
146  nn(mpi),
147  beta(mpi, CON, mp.get_bvect_spher()),
148  lnq(mpi),
149  gamma(mp.flat_met_spher()){
150 
151 
152  // Check of the EOS
153 // const Eos_poly_newt* p_eos_poly_newt =
154 // dynamic_cast<const Eos_poly_newt*>( &eos ) ;
155 //
156 // const Eos_incomp_newt* p_eos_incomp_newt =
157 // dynamic_cast<const Eos_incomp_newt*>( &eos ) ;
158 //
159 //
160 //
161 // if (p_eos_poly_newt != 0x0) {
162 // cout <<
163 // "Star::Star : the EOS Eos_poly_newt must not be employed"
164 // << " for a relativistic star ! " << endl ;
165 // cout << "(Use Eos_poly instead)" << endl ;
166 // abort() ;
167 // }
168 // if (p_eos_incomp_newt != 0x0) {
169 // cout <<
170 // "Star::Star : the EOS Eos_incomp_newt must not be employed"
171 // << " for a relativistic star ! " << endl ;
172 // cout << "(Use Eos_incomp instead)" << endl ;
173 // abort() ;
174 // }
175 //
176  // Pointers of derived quantities initialized to zero :
177  set_der_0x0() ;
178 
179  // All the matter quantities are initialized to zero :
180  nbar = 0 ;
181  ener = 0 ;
182  press = 0 ;
183  ent = 0 ;
184  temp = 0 ;
185  entropy = 0 ;
186  ener_euler = 0 ;
187  s_euler = 0 ;
188  gam_euler = 1. ;
192 
193  // The metric is initialized to the flat one :
194  Metric flat(mp.flat_met_spher()) ;
195  flat.cov() ;
196  gamma = flat ;
197 
198  logn = 0 ;
199  nn = 1. ;
200  nn.std_spectral_base() ;
201  beta.set_etat_zero() ;
202  lnq = 0 ;
203 
204 }
205 
206 // Copy constructor
207 // ----------------
209  : mp(et.mp),
210  nzet(et.nzet),
211  heos(et.heos),
212  ent(et.ent),
213  nbar(et.nbar),
214  ener(et.ener),
215  press(et.press),
216  temp(et.temp),
217  entropy(et.entropy),
218  ener_euler(et.ener_euler),
219  s_euler(et.s_euler),
220  gam_euler(et.gam_euler),
221  u_euler(et.u_euler),
222  stress_euler(et.stress_euler),
223  logn(et.logn),
224  nn(et.nn),
225  beta(et.beta),
226  lnq(et.lnq),
227  gamma(et.gamma){
228 
229  set_der_0x0() ;
230 
231 }
232 
233 // Constructor from a file
234 // -----------------------
235 Hot_star::Hot_star(Map& mpi, const Hot_eos& heos_i, FILE* fich)
236  : mp(mpi),
237  heos(heos_i),
238  ent(mpi),
239  nbar(mpi),
240  ener(mpi),
241  press(mpi),
242  temp(mpi),
243  entropy(mpi),
244  ener_euler(mpi),
245  s_euler(mpi),
246  gam_euler(mpi),
247  u_euler(mpi, CON, mp.get_bvect_spher()),
248  stress_euler(mpi, CON, mp.get_bvect_spher()),
249  logn(mpi, *(mpi.get_mg()), fich),
250  nn(mpi),
251  beta(mpi, CON, mp.get_bvect_spher()),
252  lnq(mpi, *(mpi.get_mg()), fich),
253  gamma(mpi.flat_met_spher()){
254 
255  // Hot_star parameters
256  // -----------------
257 
258  // nzet is read in the file:
259  int xx ;
260  fread_be(&xx, sizeof(int), 1, fich) ;
261  nzet = xx ;
262 
263  // Equation of state
264  // -----------------
265 
266  // Read of the saved EOS
267  Hot_eos* p_heos_file = Hot_eos::hoteos_from_file(fich) ;
268 
269  // Comparison with the assigned EOS:
270  if (heos != *p_heos_file) {
271  cout <<
272  "Hot_star::Hot_star(const Map&, const Hot_eos&, FILE*) : the hot EOS given in "
273  << endl <<
274  " argument and that read in the file are different !" << endl ;
275  abort() ;
276  }
277 
278  // p_heos_file is no longer required (it was used only for checking the
279  // EOS compatibility)
280  delete p_heos_file ;
281 
282  // Read of the saved fields:
283  // ------------------------
284  Scalar ent_file(mp, *(mp.get_mg()), fich) ;
285  ent = ent_file ;
286  Scalar entropy_file(mp, *(mp.get_mg()), fich) ;
287  entropy = entropy_file ;
290  nn = 1 ;
291  beta.set_etat_zero() ;
292 
293  // Pointers of derived quantities initialized to zero
294  // --------------------------------------------------
295  set_der_0x0() ;
296 
297 }
298 
299  //------------//
300  // Destructor //
301  //------------//
302 
304 
305  del_deriv() ;
306 
307 }
308 
309 
310  //----------------------------------//
311  // Management of derived quantities //
312  //----------------------------------//
313 
314 void Hot_star::del_deriv() const {
315 
316  if (p_mass_b != 0x0) delete p_mass_b ;
317  if (p_mass_g != 0x0) delete p_mass_g ;
318  if (p_ray_eq != 0x0) delete p_ray_eq ;
319  if (p_ray_eq_pis2 != 0x0) delete p_ray_eq_pis2 ;
320  if (p_ray_eq_pi != 0x0) delete p_ray_eq_pi ;
321  if (p_ray_eq_3pis2 != 0x0) delete p_ray_eq_3pis2 ;
322  if (p_ray_pole != 0x0) delete p_ray_pole ;
323  if (p_l_surf != 0x0) delete p_l_surf ;
324  if (p_xi_surf != 0x0) delete p_xi_surf ;
325 
327 }
328 
329 
330 
331 
332 void Hot_star::set_der_0x0() const {
333 
334  p_mass_b = 0x0 ;
335  p_mass_g = 0x0 ;
336  p_ray_eq = 0x0 ;
337  p_ray_eq_pis2 = 0x0 ;
338  p_ray_eq_pi = 0x0 ;
339  p_ray_eq_3pis2 = 0x0 ;
340  p_ray_pole = 0x0 ;
341  p_l_surf = 0x0 ;
342  p_xi_surf = 0x0 ;
343 
344 }
345 
347 
353 
354  del_deriv() ;
355 
356 }
357 
358 
359 
360 
361  //--------------//
362  // Assignment //
363  //--------------//
364 
365 // Assignment to another Hot_star
366 // ----------------------------
367 void Hot_star::operator=(const Hot_star& et) {
368 
369  assert( &(et.mp) == &mp ) ; // Same mapping
370  assert( &(et.heos) == &heos ) ; // Same EOS
371 
372  nzet = et.nzet ;
373  ent = et.ent ;
374  nbar = et.nbar ;
375  ener = et.ener ;
376  press = et.press ;
377  temp = et.temp ;
378  entropy = et.entropy ;
379  ener_euler = et.ener_euler ;
380  s_euler = et.s_euler ;
381  gam_euler = et.gam_euler ;
382  u_euler = et.u_euler ;
384  logn = et.logn ;
385  nn = et.nn ;
386  beta = et.beta ;
387  lnq = et.lnq ;
388  gamma = et.gamma ;
389 
390  del_deriv() ; // Deletes all derived quantities
391 
392 }
393 
394 // Assignment of the enthalpy field
395 // --------------------------------
396 
397 void Hot_star::set_enthalpy(const Scalar& ent_i) {
398 
399  ent = ent_i ;
400 
401  // Update of (nbar, ener, press) :
402  equation_of_state() ;
403 
404  // The derived quantities are obsolete:
405  del_deriv() ;
406 
407 }
408 
409  //--------------//
410  // Outputs //
411  //--------------//
412 
413 // Save in a file
414 // --------------
415 void Hot_star::sauve(FILE* fich) const {
416 
417  logn.sauve(fich) ;
418  lnq.sauve(fich) ;
419 
420  int xx = nzet ;
421  fwrite_be(&xx, sizeof(int), 1, fich) ;
422 
423  heos.sauve(fich) ;
424  ent.sauve(fich) ;
425  entropy.sauve(fich) ;
426 }
427 
428 // Printing
429 // --------
430 
431 ostream& operator<<(ostream& ost, const Hot_star& et) {
432  et >> ost ;
433  return ost ;
434 }
435 
436 ostream& Hot_star::operator>>(ostream& ost) const {
437 
438  using namespace Unites ;
439 
440  ost << endl ;
441 
442  ost << "Number of domains occupied by the star : " << nzet << endl ;
443 
444  ost << "Equation of state : " << endl ;
445  ost << heos << endl ;
446 
447  ost << endl << "Central enthalpy : " << ent.val_grid_point(0,0,0,0) << " c^2" << endl ;
448  ost << "Central proper baryon density : " << nbar.val_grid_point(0,0,0,0)
449  << " x 0.1 fm^-3" << endl ;
450  ost << "Central proper energy density : " << ener.val_grid_point(0,0,0,0)
451  << " rho_nuc c^2" << endl ;
452  ost << "Central pressure : " << press.val_grid_point(0,0,0,0)
453  << " rho_nuc c^2" << endl ;
454 
455  ost << endl ;
456  ost << "Central entropy per baryon : " << entropy.val_grid_point(0,0,0,0)
457  << " k_B" << endl ;
458  ost << "Total entropy (if constant) : " << entropy.val_grid_point(0,0,0,0)*mass_b()/msol
459  << " M_sol" << endl ;
460  ost << "Central temperature : " << 10/mevpfm3*temp.val_grid_point(0,0,0,0)
461  << " MeV" << endl ;
462 
463  ost << endl ;
464  ost << "Central lapse N : " << nn.val_grid_point(0,0,0,0) << endl ;
465 // ost << "Central value of lnq : " << lnq.val_grid_point(0,0,0,0) << endl ;
466 
467  ost << endl
468  << "Coordinate equatorial radius (phi=0) a1 = "
469  << ray_eq()/km << " km" << endl ;
470  ost << "Coordinate equatorial radius (phi=pi/2) a2 = "
471  << ray_eq_pis2()/km << " km" << endl ;
472  ost << "Coordinate equatorial radius (phi=pi): "
473  << ray_eq_pi()/km << " km" << endl ;
474  ost << "Coordinate polar radius a3 = "
475  << ray_pole()/km << " km" << endl ;
476  ost << "Axis ratio a2/a1 = " << ray_eq_pis2() / ray_eq()
477  << " a3/a1 = " << ray_pole() / ray_eq() << endl ;
478  ost << endl << "Baryon mass : " << mass_b() / msol << " M_sol" << endl ;
479  ost << "Gravitational mass : " << mass_g() / msol << " M_sol" << endl ;
480 
481 
482  return ost ;
483 }
484 
485  //-----------------------------------------//
486  // Computation of hydro quantities //
487  //-----------------------------------------//
488 
490 
491  Scalar ent_eos = ent ;
492 
493  // Slight rescale of the enthalpy field in case of 2 domains inside the
494  // star
495 
496  // TODO: not needed for entropy?
497 
498  double epsilon = 1.e-12 ;
499  const Mg3d* mg = mp.get_mg() ;
500  int nz = mg->get_nzone() ;
501 
502  Mtbl xi(mg) ;
503  xi.set_etat_qcq() ;
504  for (int l=0; l<nz; l++) {
505  xi.t[l]->set_etat_qcq() ;
506  for (int k=0; k<mg->get_np(l); k++) {
507  for (int j=0; j<mg->get_nt(l); j++) {
508  for (int i=0; i<mg->get_nr(l); i++) {
509  xi.set(l,k,j,i) =
510  mg->get_grille3d(l)->x[i] ;
511  }
512  }
513  }
514 
515  }
516 
517  Scalar fact_ent(mp) ;
518  fact_ent.allocate_all() ;
519 
520  fact_ent.set_domain(0) = 1 + epsilon * xi(0) * xi(0) ;
521  fact_ent.set_domain(1) = 1 - 0.25 * epsilon * (xi(1) - 1) * (xi(1) - 1) ;
522 
523  for (int l=nzet; l<nz; l++) {
524  fact_ent.set_domain(l) = 1 ;
525  }
526 
527  if (nzet > 1) {
528 
529  if(nzet == 3) {
530  fact_ent.set_domain(1) = 1 - 0.5 * epsilon * (xi(1) - 0.5) * (xi(1) - 0.5) ;
531  fact_ent.set_domain(2) = 1 - 0.25 * epsilon * (xi(2) - 1) * (xi(2) - 1) ;
532  }
533 
534  if (nzet > 3) {
535 
536  cout << "Hot_star::equation_of_state: not ready yet for nzet > 3 !"
537  << endl ;
538  }
539 
540  ent_eos = fact_ent * ent_eos ;
541  ent_eos.std_spectral_base() ;
542  }
543 
544 
545 
546  // Call to the hot EOS
547  nbar.set_etat_qcq() ;
548  nbar = heos.nbar_Hs(ent_eos, entropy, nzet) ;
549 
550  ener.set_etat_qcq() ;
551  ener = heos.ener_Hs(ent_eos, entropy, nzet) ;
552 
553  press.set_etat_qcq() ;
554  press = heos.press_Hs(ent_eos, entropy, nzet) ;
555 
556  temp.set_etat_qcq() ;
557  temp = heos.temp_Hs(ent_eos, entropy, nzet) ;
558 
559  // Set the bases for spectral expansion
564 
565  // The derived quantities are obsolete
566  del_deriv() ;
567 
568 }
569 
571 
572  cout <<
573  "Hot_star::hydro_euler : hydro_euler must be called via a derived class"
574  << endl << " of Hot_star !" << endl ;
575 
576  abort() ;
577 
578 }
579 }
Metric for tensor calculation.
Definition: metric.h:90
virtual double mass_b() const =0
Baryon mass.
Metric gamma
3-metric
Definition: hot_star.h:140
const Grille3d * get_grille3d(int l) const
Returns a pointer on the 3D mono-grid for domain no. l.
Definition: grilles.h:519
Scalar temp
Temperature field (in MeV)
Definition: hot_star.h:98
Scalar nbar
Baryon density in the fluid frame.
Definition: hot_star.h:93
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
double * p_ray_pole
Coordinate radius at .
Definition: hot_star.h:159
Multi-domain array.
Definition: mtbl.h:118
double * p_mass_g
Gravitational mass.
Definition: hot_star.h:174
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
Base class for hot stars.
Definition: hot_star.h:76
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:402
virtual void set_etat_nondef()
Sets the logical state of all components to ETATNONDEF (undefined state).
Definition: tensor.C:499
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition: hot_star.h:103
Base class for coordinate mappings.
Definition: map.h:696
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: hot_star.C:489
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
double ray_eq_pi() const
Coordinate radius at , [r_unit].
Scalar entropy
Entropy per baryon field (in $k_B$)
Definition: hot_star.h:99
virtual void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition: scalar.C:371
virtual void set_etat_nondef()
Sets the logical state to ETATNONDEF (undefined).
Definition: scalar.C:350
double * x
Array of values of at the nr collocation points.
Definition: grilles.h:212
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
Scalar nn
Lapse function N .
Definition: hot_star.h:130
void set_enthalpy(const Scalar &)
Assignment of the enthalpy field.
Definition: hot_star.C:397
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:630
Scalar press_Hs(const Scalar &ent, const Scalar &sb, int nzet, int l_min=0) const
Computes the pressure from the log-enthalpy and entropy per baryon.
Definition: hoteos.C:266
Scalar temp_Hs(const Scalar &ent, const Scalar &sb, int nzet, int l_min=0) const
Computes the temperature field from the log-enthalpy field and entropy per baryon.
Definition: hoteos.C:280
Scalar ener_Hs(const Scalar &ent, const Scalar &sb, int nzet, int l_min=0) const
Computes the total energy density from the log-enthalpy and entropy per baryon.
Definition: hoteos.C:252
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:652
Scalar logn
Logarithm of the lapse N .
Definition: hot_star.h:127
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
double ray_eq_pis2() const
Coordinate radius at , [r_unit].
double ray_eq() const
Coordinate radius at , [r_unit].
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
int nzet
Number of domains of *mp occupied by the star.
Definition: hot_star.h:84
virtual void sauve(FILE *) const
Save in a file.
Definition: hot_star.C:415
Itbl * p_l_surf
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surfa...
Definition: hot_star.h:165
double * p_ray_eq_pis2
Coordinate radius at , .
Definition: hot_star.h:150
Tbl * p_xi_surf
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the su...
Definition: hot_star.h:171
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
double * p_mass_b
Baryon mass.
Definition: hot_star.h:173
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: hot_star.h:109
Scalar ener
Total energy density in the fluid frame.
Definition: hot_star.h:94
Scalar ent
Log-enthalpy.
Definition: hot_star.h:91
Map & mp
Mapping associated with the star.
Definition: hot_star.h:81
Scalar nbar_Hs(const Scalar &ent, const Scalar &sb, int nzet, int l_min=0) const
Computes the baryon density field from the log-enthalpy field and entropy per baryon.
Definition: hoteos.C:236
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition: fwrite_be.C:73
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: mtbl.C:302
Tbl & set(int l)
Read/write of the Tbl in a given domain.
Definition: mtbl.h:221
void operator=(const Hot_star &)
Assignment to another Hot_star.
Definition: hot_star.C:367
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: hot_star.C:346
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:72
double * p_ray_eq_pi
Coordinate radius at , .
Definition: hot_star.h:153
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
Definition: metric.C:283
const Hot_eos & heos
Equation of state of the stellar matter.
Definition: hot_star.h:86
Base class for 2-parameters equations of state (abstract class).
Definition: hoteos.h:85
Multi-domain grid.
Definition: grilles.h:276
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: hot_star.C:570
Hot_star(Map &mp_i, int nzet_i, const Hot_eos &heos_i)
Standard constructor.
Definition: hot_star.C:130
double * p_ray_eq_3pis2
Coordinate radius at , .
Definition: hot_star.h:156
virtual void sauve(FILE *) const
Save in a file.
Definition: hoteos.C:129
static Hot_eos * hoteos_from_file(FILE *)
Construction of an EOS from a binary file.
Sym_tensor stress_euler
Spatial part of the stress-energy tensor with respect to the Eulerian observer.
Definition: hot_star.h:117
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition: hot_star.h:106
virtual double mass_g() const =0
Gravitational mass.
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:507
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: hot_star.C:332
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:476
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: hot_star.C:436
virtual ~Hot_star()
Destructor.
Definition: hot_star.C:303
double ray_pole() const
Coordinate radius at [r_unit].
Tbl ** t
Array (size nzone ) of pointers on the Tbl &#39;s.
Definition: mtbl.h:132
double * p_ray_eq
Coordinate radius at , .
Definition: hot_star.h:147
Scalar press
Fluid pressure.
Definition: hot_star.h:95
Vector beta
Shift vector.
Definition: hot_star.h:133
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:327
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: hot_star.C:314
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: hot_star.h:112