30 #include "star_rot_diff.h" 33 #include "utilitaires.h" 45 double (*omega_frot_i)(
double,
const Tbl&),
46 double (*primfrot_i)(
double,
const Tbl&),
47 const Tbl& par_frot_i)
48 :
Star_rot(mpi, nzet_i, relat, eos_i),
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)
84 double (*omega_frot_i)(
double,
const Tbl&),
85 double (*primfrot_i)(
double,
const Tbl&) )
87 omega_frot(omega_frot_i),
168 ost <<
"Differentially rotating star" << endl ;
169 ost <<
"-----------------------------" << endl ;
171 ost << endl <<
"Parameters of Omega(F) : " << endl ;
174 ost <<
"Min, Max of Omega/(2pi) : " <<
omega_min / (2*M_PI) * f_unit
175 <<
" Hz, " <<
omega_max / (2*M_PI) * f_unit <<
" Hz" << endl ;
176 int lsurf =
nzet - 1;
179 ost <<
"Central, equatorial value of Omega: " 183 ost <<
"Central, equatorial value of Omega/(2 Pi): " 191 ost <<
"Max prop. bar. dens. : " << nbar_max
194 ost <<
"Max prop. ener. dens. (e_max) : " << ener_max
197 ost <<
"Max pressure : " << press_max
202 double len_rho = 1. /
sqrt( ggrav * ener_max ) ;
203 ost << endl <<
"Value of A = par_frot(1) in units of e_max^{1/2} : " <<
206 ost <<
"Value of A / r_eq : " <<
209 ost <<
"r_p/r_eq : " <<
aplat() << endl ;
210 ost <<
"KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
213 /
pow(ggrav, 1.3333333333333333) << endl ;
215 ost <<
"M e_max^{1/2} : " << ggrav *
mass_g() / len_rho << endl ;
217 ost <<
"r_eq e_max^{1/2} : " <<
ray_eq() / len_rho << endl ;
219 ost <<
"T/W : " <<
tsw() << endl ;
221 ost <<
"Omega_c / e_max^{1/2} : " <<
get_omega_c() * len_rho << endl ;
240 if (p_eos_poly != 0x0) {
242 double kappa = p_eos_poly->
get_kap() ;
243 double gam = p_eos_poly->
get_gam() ; ;
246 double kap_ns2 =
pow( kappa, 0.5 /(gam-1) ) ;
249 double r_poly = kap_ns2 /
sqrt(ggrav) ;
252 double t_poly = r_poly ;
255 double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
258 ost <<
" n_max : " <<
max(
max(
nbar ) ) / rho_poly << endl ;
259 ost <<
" e_max : " <<
max(
max(
ener ) ) / rho_poly << endl ;
260 ost <<
" P_min : " << 2.*M_PI /
omega_max / t_poly << endl ;
261 ost <<
" P_max : " << 2.*M_PI /
omega_min / t_poly << endl ;
263 int lsurf =
nzet - 1;
266 ost <<
" P_eq : " << 2.*M_PI /
virtual double mass_g() const
Gravitational mass.
virtual ~Star_rot_diff()
Destructor.
void operator=(const Star_rot_diff &)
Assignment to another Star_rot_diff.
Class for differentially rotating stars in quasi-isotropic gauge and maximal slicing.
Map & mp
Mapping associated with the star.
Cmp sqrt(const Cmp &)
Square root.
Tbl par_frot
Parameters of the function .
const Eos & eos
Equation of state of the stellar matter.
Standard units of space, time and mass.
Equation of state base class.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Tensor field of valence 0 (or component of a tensorial field).
Base class for coordinate mappings.
void fait_prim_field()
Computes the member prim_field from omega_field .
double omega_funct(double F) const
Evaluates , where F is linked to the components of the fluid 4-velocity by .
double omega
Rotation angular velocity ([f_unit] )
Scalar nbar
Baryon density in the fluid frame.
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Class for isolated rotating stars.
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
void operator=(const Star_rot &)
Assignment to another Star_rot.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
int nzet
Number of domains of *mp occupied by the star.
virtual void sauve(FILE *) const
Save in a file.
double prim_funct_omega(double F) const
Evaluates the primitive of .
virtual void sauve(FILE *) const
Save in a file.
Scalar ener
Total energy density in the fluid frame.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Scalar press
Fluid pressure.
virtual double mass_b() const
Baryon mass.
double get_kap() const
Returns the pressure coefficient (cf.
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
virtual double angu_mom() const
Angular momentum.
Polytropic equation of state (relativistic case).
virtual void display_poly(ostream &) const
Display in polytropic units.
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
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.
Cmp pow(const Cmp &, int)
Power .
virtual double tsw() const
Ratio T/W.
double omega_max
Maximum value of .
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.
Star_rot_diff(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, double(*omega_frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
virtual void display_poly(ostream &) const
Display in polytropic units.
double(* omega_frot)(double, const Tbl &)
Function defining the rotation profile.
double ray_eq() const
Coordinate radius at , [r_unit].
void sauve(FILE *) const
Save in a file.
virtual void sauve(FILE *) const
Save in a file.
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
double omega_min
Minimum value of .
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
virtual double aplat() const
Flatening r_pole/r_eq.
Scalar omega_field
Field .