29 #include "star_rot_diff_cfc.h" 32 #include "utilitaires.h" 44 double (*omega_frot_i)(
double,
const Tbl&),
45 double (*primfrot_i)(
double,
const Tbl&),
46 const Tbl& par_frot_i)
48 omega_frot(omega_frot_i),
70 omega_frot(star.omega_frot),
71 primfrot(star.primfrot),
72 par_frot(star.par_frot),
73 omega_field(star.omega_field),
74 omega_min(star.omega_min),
75 omega_max(star.omega_max),
76 prim_field(star.prim_field)
82 double (*omega_frot_i)(
double,
const Tbl&),
83 double (*primfrot_i)(
double,
const Tbl&) )
85 omega_frot(omega_frot_i),
106 double lambda = 1./3. ;
178 ost <<
"Differentially rotating star" << endl ;
179 ost <<
"-----------------------------" << endl ;
181 ost << endl <<
"Parameters of Omega(F) : " << endl ;
184 ost <<
"Min, Max of Omega/(2pi) : " <<
omega_min / (2*M_PI) * f_unit
185 <<
" Hz, " <<
omega_max / (2*M_PI) * f_unit <<
" Hz" << endl ;
186 int lsurf =
nzet - 1;
189 ost <<
"Central, equatorial value of Omega: " 193 ost <<
"Central, equatorial value of Omega/(2 Pi): " 201 ost <<
"Max prop. bar. dens. : " << nbar_max
204 ost <<
"Max prop. ener. dens. (e_max) : " << ener_max
207 ost <<
"Max pressure : " << press_max
212 double len_rho = 1. /
sqrt( ggrav * ener_max ) ;
213 ost << endl <<
"Value of A = par_frot(1) in units of e_max^{1/2} : " <<
216 ost <<
"Value of A / r_eq : " <<
219 ost <<
"r_p/r_eq : " <<
aplat() << endl ;
220 ost <<
"KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
223 /
pow(ggrav, 1.3333333333333333) << endl ;
225 ost <<
"M e_max^{1/2} : " << ggrav *
mass_g() / len_rho << endl ;
227 ost <<
"r_eq e_max^{1/2} : " <<
ray_eq() / len_rho << endl ;
229 ost <<
"T/W : " <<
tsw() << endl ;
231 ost <<
"Omega_c / e_max^{1/2} : " <<
get_omega_c() * len_rho << endl ;
242 void Star_rot_diff_CFC::display_poly(ostream& ost)
const {
248 if (p_eos_poly != 0x0) {
250 double kappa = p_eos_poly->
get_kap() ;
253 double kap_ns2 =
pow( kappa, 0.5 /(p_eos_poly->
get_gam()-1) ) ;
256 double r_poly = kap_ns2 /
sqrt(ggrav) ;
259 double t_poly = r_poly ;
262 double m_poly = r_poly / ggrav ;
265 double j_poly = r_poly * r_poly / ggrav ;
268 double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
271 ost << endl <<
"Quantities in polytropic units : " << endl ;
272 ost <<
"==============================" << endl ;
273 ost <<
" ( r_poly = " << r_poly / km <<
" km )" << endl ;
276 ost <<
" Omega_c : " <<
get_omega_c() * t_poly << endl ;
277 ost <<
" P_c : " << 2.*M_PI /
get_omega_c() / t_poly << endl ;
278 ost <<
" M_bar : " <<
mass_b() / m_poly << endl ;
279 ost <<
" M : " <<
mass_g() / m_poly << endl ;
280 ost <<
" J : " <<
angu_mom() / j_poly << endl ;
281 ost <<
" r_eq : " <<
ray_eq() / r_poly << endl ;
282 ost <<
" R_circ_eq : " <<
r_circ() / r_poly << endl ;
284 ost <<
" n_max : " <<
max(
max(
nbar ) ) / rho_poly << endl ;
285 ost <<
" e_max : " <<
max(
max(
ener ) ) / rho_poly << endl ;
286 ost <<
" P_min : " << 2.*M_PI /
omega_max / t_poly << endl ;
287 ost <<
" P_max : " << 2.*M_PI /
omega_min / t_poly << endl ;
289 int lsurf =
nzet - 1;
292 ost <<
" P_eq : " << 2.*M_PI /
327 double freq = omega_c / (2.*M_PI) ;
328 ost <<
"Central Omega : " << omega_c * f_unit
329 <<
" rad/s f : " << freq * f_unit <<
" Hz" << endl ;
330 ost <<
"Rotation period : " << 1000. / (freq * f_unit) <<
" ms" 332 ost << endl <<
"Central enthalpy : " <<
ent.
val_grid_point(0,0,0,0) <<
" c^2" << endl ;
334 <<
" x 0.1 fm^-3" << endl ;
336 <<
" rho_nuc c^2" << endl ;
338 <<
" rho_nuc c^2" << endl ;
340 ost <<
"Central value of log(N) : " 345 double nphi_c =
beta(3).val_grid_point(0, 0, 0, 0) ;
346 if ( (omega_c==0) && (nphi_c==0) ) {
347 ost <<
"Central azimuthal shift : " << nphi_c << endl ;
350 ost <<
"Central azimuthal shift /Omega : " << nphi_c / omega_c
355 int lsurf =
nzet - 1;
358 ost <<
"Equatorial value of the velocity U: " 359 <<
u_euler(3).val_grid_point(
nzet-1, 0, nt-1, nr-1) <<
" c" << endl ;
362 <<
"Coordinate equatorial radius r_eq = " 363 <<
ray_eq()/km <<
" km" << endl ;
364 ost <<
"Flattening r_pole/r_eq : " <<
aplat() << endl ;
void operator=(const Star_rot_diff_CFC &)
Assignment to another Star_rot_diff.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Map & mp
Mapping associated with the star.
Cmp sqrt(const Cmp &)
Square root.
double omega_min
Minimum value of .
Scalar psi
Conformal factor .
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).
Scalar omega_field
Field .
virtual double angu_mom() const
Angular momentum.
Base class for coordinate mappings.
virtual double r_circ() const
Circumferential equatorial radius.
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
void operator=(const Star_rot_CFC &)
Assignment to another Star_rot_CFC.
Scalar psi4
Conformal factor .
Class for differentially rotating stars in Conformal Flatness Condition and maximal slicing...
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Scalar nbar
Baryon density in the fluid frame.
virtual void sauve(FILE *) const
Save in a file.
Tensor field of valence 1.
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
virtual ~Star_rot_diff_CFC()
Destructor.
int nzet
Number of domains of *mp occupied by the star.
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
Class for relativistic rotating stars in Conformal Flatness Condition and maximal slicing...
virtual void sauve(FILE *) const
Save in a file.
double prim_funct_omega(double F) const
Evaluates the primitive of .
Scalar ener
Total energy density in the fluid frame.
virtual double aplat() const
Flattening r_pole/r_eq.
double omega
Rotation angular velocity ([f_unit] )
Scalar press
Fluid pressure.
double get_kap() const
Returns the pressure coefficient (cf.
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
virtual void sauve(FILE *) const
Save in a file.
Polytropic equation of state (relativistic case).
Tbl par_frot
Parameters of the function .
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
virtual double mass_b() const
Baryonic mass.
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
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 .
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Scalar logn
Logarithm of the lapse N .
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.
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
double(* omega_frot)(double, const Tbl &)
Function defining the rotation profile.
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Scalar nn
Lapse function N .
double omega_max
Maximum value of .
Star_rot_diff_CFC(Map &mp_i, int nzet_i, const Eos &eos_i, int relat_i, int filter, double(*omega_frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
double ray_eq() const
Coordinate radius at , [r_unit].
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
void sauve(FILE *) const
Save in a file.
virtual double tsw() const
Ratio T/W.
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
double omega_funct(double F) const
Evaluates , where F is linked to the components of the fluid 4-velocity by .
Vector j_euler
Momentum density 3-vector with respect to the Eulerian observer.
void fait_prim_field()
Computes the member prim_field from omega_field .
virtual double mass_g() const
Gravitational mass.
const Metric_flat & flat
flat metric (spherical components)