|
LORENE
|
Class for computing a Black hole - Neutron star binary system with an extreme mass ratio. More...
#include <bin_bhns_extr.h>
Public Member Functions | |
| Bin_bhns_extr (Map &mp, int nzet, const Eos &eos, bool irrot, bool relat, bool kerrs, bool multi) | |
| Standard constructor. | |
| Bin_bhns_extr (const Bin_bhns_extr &) | |
| Copy constructor. | |
| Bin_bhns_extr (Map &mp, const Eos &eos, FILE *fich) | |
Constructor from a file (see sauve(FILE*) ). | |
| ~Bin_bhns_extr () | |
| Destructor. | |
| void | operator= (const Bin_bhns_extr &) |
| Assignment to another Bin_bhns_extr. | |
| Et_bin_bhns_extr & | set_ns () |
| Read/write of the neutron star. | |
| double & | set_omega () |
| Sets the orbital angular velocity [{\tt f_unit}]. | |
| double & | set_separ () |
| Sets the orbital separation [{\tt r_unit}]. | |
| double & | set_mass_bh () |
| Sets the gravitational mass of BH [{\tt m_unit}]. | |
| const Et_bin_bhns_extr & | get_ns () const |
| Returns a reference to the neutron star. | |
| double | get_omega () const |
| Returns the orbital angular velocity [{\tt f_unit}]. | |
| double | get_separ () const |
| Returns the coordinate separation of the binary system [{\tt r_unit}]. | |
| double | get_mass_bh () const |
| Returns the gravitational mass of BH [{\tt m_unit}]. | |
| void | sauve (FILE *) const |
| Save in a file. | |
| void | display_poly (ostream &) const |
| Display in polytropic units. | |
| double | xa_barycenter_extr () const |
| Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or in the conformally flat one. | |
| double | ya_barycenter_extr () const |
| in the Kerr-Schild background metric | |
| double | mass_b_extr () const |
| Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat. | |
| void | orbit_omega_ks (double fact_omeg_min, double fact_omeg_max) |
| Computes the orbital angular velocity {\tt omega} in the Kerr-Schild background metric. | |
| void | orbit_omega_cf (double fact_omeg_min, double fact_omeg_max) |
| Computes the orbital angular velocity {\tt omega} in the conformally flat background metric. | |
| void | analytical_omega () |
| Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian case). | |
| void | analytical_shift () |
| Sets some analytical template for the shift vector (via the members {\tt w_shift} and {\tt khi_shift} of {\tt Etoile_bin}). | |
Private Member Functions | |
| void | del_deriv () const |
| Deletes all the derived quantities. | |
| void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. | |
| ostream & | operator>> (ostream &) const |
| Operator >> (function called by the operator <<). | |
Private Attributes | |
| const Base_vect_cart | ref_triad |
| Cartesian triad of the absolute reference frame. | |
| Et_bin_bhns_extr | star |
| Neutron star. | |
| double | omega |
| Angular velocity with respect to an asymptotically inertial observer. | |
| double | separ |
| Absolute orbital separation between two centers of BH and NS. | |
| double | mass_bh |
| Gravitational mass of BH. | |
| double * | p_xa_barycenter_extr |
| Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or in the conformally flat one. | |
| double * | p_ya_barycenter_extr |
| Absolute coordinate Y of the barycenter of the baryon density in the Kerr-Schild background metric. | |
| double * | p_mass_b_extr |
| Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat one. | |
Friends | |
| ostream & | operator<< (ostream &ost, const Bin_bhns_extr &bibi) |
| Display. | |
Class for computing a Black hole - Neutron star binary system with an extreme mass ratio.
()
Definition at line 57 of file bin_bhns_extr.h.
| Lorene::Bin_bhns_extr::Bin_bhns_extr | ( | Map & | mp, |
| int | nzet, | ||
| const Eos & | eos, | ||
| bool | irrot, | ||
| bool | relat, | ||
| bool | kerrs, | ||
| bool | multi ) |
Standard constructor.
| mp | Mapping on which the neutron star will be defined |
| nzet | Number of domains occupied by the neutron star |
| eos | Equation of state of the neutron star |
| irrot | should be {\tt true} if NS is irrotational, {\tt false} if NS is corotating |
| relat | should be {\tt true} for a relativistic configuration, {\tt false} for a Newtonian one |
| kerrs | should be {\tt true} for the Kerr-Schild background metric, {\tt false} for the conformally flat one |
| multi | should be {\tt true} for the multipole falloff boundary condition, {\tt true} for the ![]() |
Definition at line 85 of file bin_bhns_extr.C.
References Lorene::Map(), mass_bh, omega, ref_triad, separ, set_der_0x0(), and star.
| Lorene::Bin_bhns_extr::Bin_bhns_extr | ( | const Bin_bhns_extr & | bibi | ) |
Copy constructor.
Definition at line 102 of file bin_bhns_extr.C.
References Bin_bhns_extr(), mass_bh, omega, ref_triad, separ, set_der_0x0(), and star.
Constructor from a file (see sauve(FILE*) ).
| mp | Mapping on which the neutron star will be defined |
| eos | Equation of state of the neutron star |
Definition at line 117 of file bin_bhns_extr.C.
References Lorene::fread_be(), Lorene::Map(), mass_bh, omega, ref_triad, separ, set_der_0x0(), and star.
| Lorene::Bin_bhns_extr::~Bin_bhns_extr | ( | ) |
| void Lorene::Bin_bhns_extr::analytical_omega | ( | ) |
Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian case).
Definition at line 55 of file bin_bhns_extr_omegaana.C.
References del_deriv(), mass_bh, omega, Lorene::pow(), separ, Lorene::sqrt(), and star.
| void Lorene::Bin_bhns_extr::analytical_shift | ( | ) |
Sets some analytical template for the shift vector (via the members {\tt w_shift} and {\tt khi_shift} of {\tt Etoile_bin}).
Definition at line 55 of file bin_bhns_extr_anashift.C.
References Lorene::Cmp::annule(), Lorene::Map(), omega, separ, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), and star.
|
private |
Deletes all the derived quantities.
Definition at line 145 of file bin_bhns_extr.C.
References p_mass_b_extr, p_xa_barycenter_extr, p_ya_barycenter_extr, and set_der_0x0().
| void Lorene::Bin_bhns_extr::display_poly | ( | ostream & | ost | ) | const |
Display in polytropic units.
Definition at line 323 of file bin_bhns_extr.C.
References Lorene::Eos_poly::get_gam(), Lorene::Eos_poly::get_kap(), mass_b_extr(), mass_bh, omega, Lorene::pow(), separ, Lorene::sqrt(), star, xa_barycenter_extr(), and ya_barycenter_extr().
|
inline |
Returns the gravitational mass of BH [{\tt m_unit}].
Definition at line 178 of file bin_bhns_extr.h.
References mass_bh.
|
inline |
Returns a reference to the neutron star.
Definition at line 166 of file bin_bhns_extr.h.
References star.
|
inline |
Returns the orbital angular velocity [{\tt f_unit}].
Definition at line 170 of file bin_bhns_extr.h.
References omega.
|
inline |
Returns the coordinate separation of the binary system [{\tt r_unit}].
Definition at line 175 of file bin_bhns_extr.h.
References separ.
| double Lorene::Bin_bhns_extr::mass_b_extr | ( | ) | const |
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat.
Definition at line 218 of file bin_bhns_extr_global.C.
References Lorene::Cmp::integrale(), Lorene::Map(), mass_bh, p_mass_b_extr, Lorene::pow(), separ, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::sqrt(), star, and Lorene::Cmp::std_base_scal().
| void Lorene::Bin_bhns_extr::operator= | ( | const Bin_bhns_extr & | bibi | ) |
Assignment to another Bin_bhns_extr.
Definition at line 169 of file bin_bhns_extr.C.
References Bin_bhns_extr(), del_deriv(), mass_bh, omega, ref_triad, separ, and star.
|
private |
Operator >> (function called by the operator <<).
Definition at line 210 of file bin_bhns_extr.C.
References mass_b_extr(), mass_bh, omega, separ, and star.
| void Lorene::Bin_bhns_extr::orbit_omega_cf | ( | double | fact_omeg_min, |
| double | fact_omeg_max ) |
Computes the orbital angular velocity {\tt omega} in the conformally flat background metric.
| fact_omeg_min | [input] : determines the lower bound of the interval {\tt [omega_min, omega_max]} in which {\tt omega} is searched by {\tt omega_min = fact_omeg_min * omega}, where {\tt omega} is the previous value of the angular velocity (typical value : {\tt fact_omeg_min = 0.5}) |
| fact_omeg_max | [input] : determines the higher bound of the interval {\tt [omega_min, omega_max]} in which {\tt omega} is searched by {\tt omega_max = fact_omeg_max * omega}, where {\tt omega} is the previous value of the angular velocity. (typical value : {\tt fact_omeg_max = 1.5}) |
Definition at line 322 of file bin_bhns_extr_orbit.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Tenseur::change_triad(), Lorene::Cmp::dsdx(), Lorene::flat_scalar_prod(), Lorene::Tbl::get_taille(), Lorene::Tenseur::get_triad(), Lorene::Map(), omega, ref_triad, star, Lorene::zero_list(), and Lorene::zerosec_b().
| void Lorene::Bin_bhns_extr::orbit_omega_ks | ( | double | fact_omeg_min, |
| double | fact_omeg_max ) |
Computes the orbital angular velocity {\tt omega} in the Kerr-Schild background metric.
| fact_omeg_min | [input] : determines the lower bound of the interval {\tt [omega_min, omega_max]} in which {\tt omega} is searched by {\tt omega_min = fact_omeg_min * omega}, where {\tt omega} is the previous value of the angular velocity (typical value : {\tt fact_omeg_min = 0.5}) |
| fact_omeg_max | [input] : determines the higher bound of the interval {\tt [omega_min, omega_max]} in which {\tt omega} is searched by {\tt omega_max = fact_omeg_max * omega}, where {\tt omega} is the previous value of the angular velocity. (typical value : {\tt fact_omeg_max = 1.5}) |
Definition at line 70 of file bin_bhns_extr_orbit.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Tenseur::change_triad(), Lorene::Cmp::dsdx(), Lorene::flat_scalar_prod(), Lorene::Tbl::get_taille(), Lorene::Tenseur::get_triad(), Lorene::Map(), mass_bh, omega, ref_triad, separ, star, Lorene::zero_list(), and Lorene::zerosec_b().
| void Lorene::Bin_bhns_extr::sauve | ( | FILE * | fich | ) | const |
Save in a file.
Definition at line 191 of file bin_bhns_extr.C.
References Lorene::fwrite_be(), mass_bh, omega, separ, and star.
|
private |
Sets to 0x0 all the pointers on derived quantities.
Definition at line 155 of file bin_bhns_extr.C.
References p_mass_b_extr, p_xa_barycenter_extr, and p_ya_barycenter_extr.
|
inline |
Sets the gravitational mass of BH [{\tt m_unit}].
Definition at line 160 of file bin_bhns_extr.h.
References mass_bh.
|
inline |
Read/write of the neutron star.
Definition at line 149 of file bin_bhns_extr.h.
References del_deriv(), and star.
|
inline |
Sets the orbital angular velocity [{\tt f_unit}].
Definition at line 154 of file bin_bhns_extr.h.
References omega.
|
inline |
Sets the orbital separation [{\tt r_unit}].
Definition at line 157 of file bin_bhns_extr.h.
References separ.
| double Lorene::Bin_bhns_extr::xa_barycenter_extr | ( | ) | const |
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or in the conformally flat one.
Definition at line 63 of file bin_bhns_extr_global.C.
References Lorene::Cmp::integrale(), Lorene::Map(), mass_b_extr(), mass_bh, p_xa_barycenter_extr, Lorene::pow(), separ, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::sqrt(), star, and Lorene::Cmp::std_base_scal().
| double Lorene::Bin_bhns_extr::ya_barycenter_extr | ( | ) | const |
in the Kerr-Schild background metric
Definition at line 140 of file bin_bhns_extr_global.C.
References Lorene::Cmp::integrale(), Lorene::Map(), mass_b_extr(), mass_bh, p_ya_barycenter_extr, Lorene::pow(), separ, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::sqrt(), star, and Lorene::Cmp::std_base_scal().
|
friend |
Display.
Definition at line 203 of file bin_bhns_extr.C.
References Bin_bhns_extr(), and operator<<.
|
private |
Gravitational mass of BH.
Definition at line 77 of file bin_bhns_extr.h.
|
private |
Angular velocity with respect to an asymptotically inertial observer.
Definition at line 71 of file bin_bhns_extr.h.
|
mutableprivate |
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat one.
Definition at line 96 of file bin_bhns_extr.h.
|
mutableprivate |
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or in the conformally flat one.
Definition at line 86 of file bin_bhns_extr.h.
|
mutableprivate |
Absolute coordinate Y of the barycenter of the baryon density in the Kerr-Schild background metric.
Definition at line 91 of file bin_bhns_extr.h.
|
private |
Cartesian triad of the absolute reference frame.
Definition at line 63 of file bin_bhns_extr.h.
|
private |
Absolute orbital separation between two centers of BH and NS.
Definition at line 74 of file bin_bhns_extr.h.
|
private |
Neutron star.
Definition at line 66 of file bin_bhns_extr.h.