29 #include "hot_star_rot_diff_cfc.h" 32 #include "utilitaires.h" 44 int e_type_i,
double const_value_i,
int filter,
45 double (*omega_frot_i)(
double,
const Tbl&),
46 double (*primfrot_i)(
double,
const Tbl&),
47 const Tbl& par_frot_i)
48 :
Hot_star_rot_CFC(mpi, nzet_i, heos_i, relat_i, e_type_i, const_value_i, filter),
49 omega_frot(omega_frot_i),
71 omega_frot(star.omega_frot),
72 primfrot(star.primfrot),
73 par_frot(star.par_frot),
74 omega_field(star.omega_field),
75 omega_min(star.omega_min),
76 omega_max(star.omega_max),
77 prim_field(star.prim_field)
83 double (*omega_frot_i)(
double,
const Tbl&),
84 double (*primfrot_i)(
double,
const Tbl&) )
86 omega_frot(omega_frot_i),
107 double lambda = 1./3. ;
179 ost <<
"Differentially rotating star" << endl ;
180 ost <<
"-----------------------------" << endl ;
182 ost << endl <<
"Parameters of Omega(F) : " << endl ;
185 ost <<
"Min, Max of Omega/(2pi) : " <<
omega_min / (2*M_PI) * f_unit
186 <<
" Hz, " <<
omega_max / (2*M_PI) * f_unit <<
" Hz" << endl ;
187 int lsurf =
nzet - 1;
190 ost <<
"Central, equatorial value of Omega: " 194 ost <<
"Central, equatorial value of Omega/(2 Pi): " 202 ost <<
"Max prop. bar. dens. : " << nbar_max
205 ost <<
"Max prop. ener. dens. (e_max) : " << ener_max
208 ost <<
"Max pressure : " << press_max
213 double len_rho = 1. /
sqrt( ggrav * ener_max ) ;
214 ost << endl <<
"Value of A = par_frot(1) in units of e_max^{1/2} : " <<
217 ost <<
"Value of A / r_eq : " <<
220 ost <<
"r_p/r_eq : " <<
aplat() << endl ;
221 ost <<
"KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
224 /
pow(ggrav, 1.3333333333333333) << endl ;
226 ost <<
"M e_max^{1/2} : " << ggrav *
mass_g() / len_rho << endl ;
228 ost <<
"r_eq e_max^{1/2} : " <<
ray_eq() / len_rho << endl ;
230 ost <<
"T/W : " <<
tsw() << endl ;
232 ost <<
"Omega_c / e_max^{1/2} : " <<
get_omega_c() * len_rho << endl ;
266 double freq = omega_c / (2.*M_PI) ;
267 ost <<
"Central Omega : " << omega_c * f_unit
268 <<
" rad/s f : " << freq * f_unit <<
" Hz" << endl ;
269 ost <<
"Rotation period : " << 1000. / (freq * f_unit) <<
" ms" 271 ost << endl <<
"Central enthalpy : " <<
ent.
val_grid_point(0,0,0,0) <<
" c^2" << endl ;
273 <<
" x 0.1 fm^-3" << endl ;
275 <<
" rho_nuc c^2" << endl ;
277 <<
" rho_nuc c^2" << endl ;
279 ost <<
"Central value of log(N) : " 284 double nphi_c =
beta(3).val_grid_point(0, 0, 0, 0) ;
285 if ( (omega_c==0) && (nphi_c==0) ) {
286 ost <<
"Central azimuthal shift : " << nphi_c << endl ;
289 ost <<
"Central azimuthal shift /Omega : " << nphi_c / omega_c
294 int lsurf =
nzet - 1;
297 ost <<
"Equatorial value of the velocity U: " 298 <<
u_euler(3).val_grid_point(
nzet-1, 0, nt-1, nr-1) <<
" c" << endl ;
301 <<
"Coordinate equatorial radius r_eq = " 302 <<
ray_eq()/km <<
" km" << endl ;
303 ost <<
"Flattening r_pole/r_eq : " <<
aplat() << endl ;
const Metric_flat & flat
flat metric (spherical components)
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
Scalar nbar
Baryon density in the fluid frame.
Scalar psi
Conformal factor .
Hot_star_rot_diff_CFC(Map &mp_i, int nzet_i, const Hot_eos &heos_i, int relat_i, int e_type, double const_value_i, int filter, double(*omega_frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
Cmp sqrt(const Cmp &)
Square root.
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
virtual double aplat() const
Flattening r_pole/r_eq.
Tensor field of valence 0 (or component of a tensorial field).
Tbl par_frot
Parameters of the function .
Base class for coordinate mappings.
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Class for hot differentially rotating stars in Conformal Flatness Condition and maximal slicing...
void operator=(const Hot_star_rot_diff_CFC &)
Assignment to another Hot_star_rot_diff_CFC.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Scalar nn
Lapse function N .
Tensor field of valence 1.
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Scalar logn
Logarithm of the lapse N .
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
double ray_eq() const
Coordinate radius at , [r_unit].
virtual void sauve(FILE *) const
Save in a file.
int nzet
Number of domains of *mp occupied by the star.
virtual void sauve(FILE *) const
Save in a file.
virtual double mass_b() const
Baryonic mass.
double omega
Rotation angular velocity ([f_unit] )
void operator=(const Hot_star_rot_CFC &)
Assignment to another Hot_star_rot_CFC.
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
double omega_funct(double F) const
Evaluates , where F is linked to the components of the fluid 4-velocity by .
virtual double mass_g() const
Gravitational mass.
Scalar ener
Total energy density in the fluid frame.
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Map & mp
Mapping associated with the star.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Scalar psi4
Conformal factor .
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.
double(* omega_frot)(double, const Tbl &)
Function defining the rotation profile.
Cmp pow(const Cmp &, int)
Power .
virtual double angu_mom() const
Angular momentum.
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Class for relativistic rotating stars in Conformal Flatness Condition and maximal slicing...
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Base class for 2-parameters equations of state (abstract class).
double omega_max
Maximum value of .
virtual void sauve(FILE *) const
Save in a file.
virtual double tsw() const
Ratio T/W.
void fait_prim_field()
Computes the member prim_field from omega_field .
void sauve(FILE *) const
Save in a file.
Scalar omega_field
Field .
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Scalar press
Fluid pressure.
Vector j_euler
Momentum density 3-vector with respect to the Eulerian observer.
double omega_min
Minimum value of .
virtual ~Hot_star_rot_diff_CFC()
Destructor.
double prim_funct_omega(double F) const
Evaluates the primitive of .
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.