LORENE
Lorene::Bin_bhns Class Reference

Class for computing a black hole - neutron star binary system with comparable mass () More...

#include <bin_bhns.h>

Public Member Functions

 Bin_bhns (Map &mp_bh, Map &mp_ns, int nzet, const Eos &eos, bool irrot_ns, bool kerrschild, bool bc_lapse_nd, bool bc_lapse_fs, bool irrot_bh, double mass_bh)
 Relative error on the Hamiltonian constraint. More...
 
 Bin_bhns (const Bin_bhns &)
 Copy constructor. More...
 
 Bin_bhns (Map &mp_bh, Map &mp_ns, const Eos &eos, FILE *fich)
 Constructor from a file (see sauve(FILE*) ) More...
 
virtual ~Bin_bhns ()
 Destructor. More...
 
void operator= (const Bin_bhns &)
 Assignment to another Bin_bhns. More...
 
Hole_bhnsset_bh ()
 Read/write of the black hole. More...
 
Star_bhnsset_ns ()
 Read/write of the neutron star. More...
 
double & set_omega ()
 Sets the orbital angular velocity [{ f_unit}]. More...
 
double & set_separ ()
 Sets the orbital separation [{ r_unit}]. More...
 
double & set_x_rot ()
 Sets the absolute coordinate X of the rotation axis [{ r_unit}]. More...
 
double & set_y_rot ()
 Sets the absolute coordinate Y of the rotation axis [{ r_unit}]. More...
 
const Hole_bhnsget_bh () const
 Returns a reference to the black hole. More...
 
const Star_bhnsget_ns () const
 Returns a reference to the neutron star. More...
 
double get_omega () const
 Returns the orbital angular velocity [{ f_unit}]. More...
 
double get_separ () const
 Returns the coordinate separation of the binary system [{ r_unit}]. More...
 
double get_x_rot () const
 Returns the absolute coordinate X of the rotation axis [{ r_unit}]. More...
 
double get_y_rot () const
 Returns the absolute coordinate Y of the rotation axis [{ r_unit}]. More...
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
void display_poly (ostream &) const
 Display in polytropic units. More...
 
double mass_adm_bhns_surf () const
 Total ADM mass. More...
 
double mass_adm_bhns_vol () const
 
double mass_kom_bhns_surf () const
 Total Komar mass. More...
 
double mass_kom_bhns_vol () const
 
const Tblline_mom_bhns () const
 Total linear momentum. More...
 
const Tblangu_mom_bhns () const
 Total angular momentum. More...
 
double virial_bhns_surf () const
 Estimates the relative error on the virial theorem $|1 - M_{ Komar} / M_{ ADM}|$. More...
 
double virial_bhns_vol () const
 Estimates the relative error on the virial theorem $|1 - M_{ Komar} / M_{ ADM}|$. More...
 
double xa_barycenter () const
 Absolute coordinate X of the barycenter of the baryon density. More...
 
double ya_barycenter () const
 Absolute coordinate Y of the barycenter of the baryon density. More...
 
double omega_two_points () const
 Orbital angular velocity derived from another method. More...
 
void orbit_omega (double fact_omeg_min, double fact_omeg_max)
 Estimates the relative error on the Hamiltonian constraint. More...
 
void rotation_axis_x (double rot_exp_x)
 Computes the position of the rotation axis X. More...
 
void rotation_axis_y (double thres_rot, double rot_exp_y, double fact)
 Computes the position of the rotation axis Y. More...
 
void shift_analytic (double reduce_shift_bh, double reduce_shift_ns)
 Sets some analytical template for the initial shift vector. More...
 

Protected Member Functions

void del_deriv () const
 Deletes all the derived quantities. More...
 
void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities. More...
 

Protected Attributes

const Base_vect_cart ref_triad
 Cartesian triad of the absolute reference frame. More...
 
Hole_bhns hole
 Black hole. More...
 
Star_bhns star
 Neutron star. More...
 
double omega
 Angular velocity with respect to an asymptotically inertial observer. More...
 
double separ
 Absolute orbital separation between two centers of BH and NS. More...
 
double x_rot
 Absolute X coordinate of the rotation axis. More...
 
double y_rot
 Absolute Y coordinate of the rotation axis. More...
 
double * p_mass_adm_bhns_surf
 Total ADM mass of the system calculated by the surface integral at infinity. More...
 
double * p_mass_adm_bhns_vol
 Total ADM mass of the system calculated by the volume integral and the surface integral at the apparent horizon. More...
 
double * p_mass_kom_bhns_surf
 Total Komar mass of the system calculated by the surface integral at infinity. More...
 
double * p_mass_kom_bhns_vol
 Total Komar mass of the system calculated by the volume integral and the surface integral at the apparent horizon. More...
 
Tblp_line_mom_bhns
 Total linear momentum of the system. More...
 
Tblp_angu_mom_bhns
 Total angular momentum of the system. More...
 
double * p_virial_bhns_surf
 Virial theorem error calculated by the ADM mass and the Komar mass of the surface integral at infinity. More...
 
double * p_virial_bhns_vol
 Virial theorem error calculated by the ADM mass and the Komar mass of the volume integral. More...
 
double * p_xa_barycenter
 Absolute coordinate X of the barycenter of the baryon density. More...
 
double * p_ya_barycenter
 Absolute coordinate Y of the barycenter of the baryon density. More...
 
double * p_omega_two_points
 Orbital angular velocity derived from another method. More...
 

Private Member Functions

ostream & operator>> (ostream &) const
 Operator >> (function called by the operator <<) More...
 

Friends

ostream & operator<< (ostream &, const Bin_bhns &)
 Display. More...
 

Detailed Description

Class for computing a black hole - neutron star binary system with comparable mass ()

Definition at line 63 of file bin_bhns.h.

Constructor & Destructor Documentation

◆ Bin_bhns() [1/3]

Lorene::Bin_bhns::Bin_bhns ( Map mp_bh,
Map mp_ns,
int  nzet,
const Eos eos,
bool  irrot_ns,
bool  kerrschild,
bool  bc_lapse_nd,
bool  bc_lapse_fs,
bool  irrot_bh,
double  mass_bh 
)

Relative error on the Hamiltonian constraint.

Relative error on the momentum constraint Standard constructor

Parameters
mp_bhMapping on which the black hole will be defined
mp_nsMapping on which the neutron star will be defined
nzetNumber of domains occupied by the neutron star
eosEquation of state of the neutron star
irrot_nsshould be { true} if NS is irrotational, { false} if NS is corotating
kerrschildshould be { true} if the background metric is Kerr-Schild, { false} if the background metric is conformally flat
bc_lapse_ndshould be { true} if the BC type for lapse is Neumann, { false} if the BC type is Dirichlet
bc_lapse_fsshould be { true} if the BC is first type { false} if the BC is second one
irrot_bhshould be { true} if BH is irrotational, { false} if NS is corotating
mass_bhBlack hole mass which appears in the background metric

Definition at line 73 of file bin_bhns.C.

References omega, separ, set_der_0x0(), x_rot, and y_rot.

◆ Bin_bhns() [2/3]

Lorene::Bin_bhns::Bin_bhns ( const Bin_bhns bibi)

Copy constructor.

Definition at line 94 of file bin_bhns.C.

References set_der_0x0().

◆ Bin_bhns() [3/3]

Lorene::Bin_bhns::Bin_bhns ( Map mp_bh,
Map mp_ns,
const Eos eos,
FILE *  fich 
)

Constructor from a file (see sauve(FILE*) )

Definition at line 111 of file bin_bhns.C.

References Lorene::fread_be(), omega, separ, set_der_0x0(), x_rot, and y_rot.

◆ ~Bin_bhns()

Lorene::Bin_bhns::~Bin_bhns ( )
virtual

Destructor.

Definition at line 133 of file bin_bhns.C.

References del_deriv().

Member Function Documentation

◆ angu_mom_bhns()

const Tbl & Lorene::Bin_bhns::angu_mom_bhns ( ) const

◆ del_deriv()

void Lorene::Bin_bhns::del_deriv ( ) const
protected

◆ display_poly()

◆ get_bh()

const Hole_bhns& Lorene::Bin_bhns::get_bh ( ) const
inline

Returns a reference to the black hole.

Definition at line 225 of file bin_bhns.h.

References hole.

◆ get_ns()

const Star_bhns& Lorene::Bin_bhns::get_ns ( ) const
inline

Returns a reference to the neutron star.

Definition at line 228 of file bin_bhns.h.

References star.

◆ get_omega()

double Lorene::Bin_bhns::get_omega ( ) const
inline

Returns the orbital angular velocity [{ f_unit}].

Definition at line 231 of file bin_bhns.h.

References omega.

◆ get_separ()

double Lorene::Bin_bhns::get_separ ( ) const
inline

Returns the coordinate separation of the binary system [{ r_unit}].

Definition at line 236 of file bin_bhns.h.

References separ.

◆ get_x_rot()

double Lorene::Bin_bhns::get_x_rot ( ) const
inline

Returns the absolute coordinate X of the rotation axis [{ r_unit}].

Definition at line 241 of file bin_bhns.h.

References x_rot.

◆ get_y_rot()

double Lorene::Bin_bhns::get_y_rot ( ) const
inline

Returns the absolute coordinate Y of the rotation axis [{ r_unit}].

Definition at line 246 of file bin_bhns.h.

References y_rot.

◆ line_mom_bhns()

◆ mass_adm_bhns_surf()

◆ mass_kom_bhns_surf()

◆ omega_two_points()

◆ operator=()

void Lorene::Bin_bhns::operator= ( const Bin_bhns bibi)

Assignment to another Bin_bhns.

Definition at line 190 of file bin_bhns.C.

References del_deriv(), hole, omega, ref_triad, separ, star, x_rot, and y_rot.

◆ operator>>()

◆ orbit_omega()

void Lorene::Bin_bhns::orbit_omega ( double  fact_omeg_min,
double  fact_omeg_max 
)

Estimates the relative error on the Hamiltonian constraint.

Estimates the relative error on the momentum constraintComputes the orbital angular velocity { omega}

fact_omeg_min [input] : determines the lower bound of the interval { [omega_min, omega_max]} in which { omega} is searched by { omega_min = fact_omeg_min * omega}, where { omega} is the previous value of the angular velocity (typical value : { fact_omeg_min = 0.5})

Parameters
fact_omeg_max[input] : determines the higher bound of the interval { [omega_min, omega_max]} in which { omega} is searched by { omega_max = fact_omeg_max * omega}, where { omega} is the previous value of the angular velocity. (typical value : { fact_omeg_max = 1.5})

Definition at line 71 of file bin_bhns_orbit.C.

References Lorene::Param::add_double(), Lorene::Scalar::dsdx(), Lorene::Star_bhns::get_confo_auto(), Lorene::Star_bhns::get_confo_tot(), Lorene::Star_bhns::get_d_confo_comp(), Lorene::Star_bhns::get_d_lapconf_comp(), Lorene::Star_bhns::get_d_shift_comp(), Lorene::Star_bhns::get_lapconf_auto(), Lorene::Star_bhns::get_lapconf_tot(), Lorene::Star_bhns::get_loggam(), Lorene::Black_hole::get_mass_bh(), Lorene::Black_hole::get_mp(), Lorene::Star::get_mp(), Lorene::Map::get_rot_phi(), Lorene::Star_bhns::get_shift_auto(), Lorene::Star_bhns::get_shift_tot(), Lorene::Tbl::get_taille(), hole, Lorene::Black_hole::is_kerrschild(), omega, Lorene::pow(), separ, star, Lorene::Scalar::std_spectral_base(), Lorene::Scalar::val_grid_point(), x_rot, y_rot, Lorene::zero_list(), and Lorene::zerosec_b().

◆ rotation_axis_x()

void Lorene::Bin_bhns::rotation_axis_x ( double  rot_exp_x)

Computes the position of the rotation axis X.

Parameters
rot_exp_x[input] : exponent of the factor which modifies the position of the two stars from the rotation axis

Definition at line 61 of file bin_bhns_rotaxis.C.

References Lorene::Black_hole::get_mass_bh(), Lorene::Black_hole::get_mp(), Lorene::Star::get_mp(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), hole, line_mom_bhns(), omega, Lorene::pow(), separ, Lorene::Black_hole::set_mp(), Lorene::Star::set_mp(), set_x_rot(), and star.

◆ rotation_axis_y()

void Lorene::Bin_bhns::rotation_axis_y ( double  thres_rot,
double  rot_exp_y,
double  fact 
)

Computes the position of the rotation axis Y.

Parameters
thres_rot[input] : threshold to stop moving to the Y dir.
rot_exp_y[input] : exponent of the factor which modifies the Y position of the neutron star coordinate
fact[input] : factor to multiply to Y_NS

Definition at line 104 of file bin_bhns_rotaxis.C.

References Lorene::Black_hole::get_mass_bh(), Lorene::Black_hole::get_mp(), Lorene::Star::get_mp(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), hole, line_mom_bhns(), omega, Lorene::pow(), separ, Lorene::Black_hole::set_mp(), Lorene::Star::set_mp(), set_y_rot(), star, xa_barycenter(), and ya_barycenter().

◆ sauve()

void Lorene::Bin_bhns::sauve ( FILE *  fich) const
virtual

Save in a file.

Definition at line 213 of file bin_bhns.C.

References Lorene::fwrite_be(), hole, omega, Lorene::Star_bhns::sauve(), Lorene::Hole_bhns::sauve(), separ, star, x_rot, and y_rot.

◆ set_bh()

Hole_bhns& Lorene::Bin_bhns::set_bh ( )
inline

Read/write of the black hole.

Definition at line 199 of file bin_bhns.h.

References del_deriv(), and hole.

◆ set_der_0x0()

void Lorene::Bin_bhns::set_der_0x0 ( ) const
protected

◆ set_ns()

Star_bhns& Lorene::Bin_bhns::set_ns ( )
inline

Read/write of the neutron star.

Definition at line 204 of file bin_bhns.h.

References del_deriv(), and star.

◆ set_omega()

double& Lorene::Bin_bhns::set_omega ( )
inline

Sets the orbital angular velocity [{ f_unit}].

Definition at line 209 of file bin_bhns.h.

References omega.

◆ set_separ()

double& Lorene::Bin_bhns::set_separ ( )
inline

Sets the orbital separation [{ r_unit}].

Definition at line 212 of file bin_bhns.h.

References separ.

◆ set_x_rot()

double& Lorene::Bin_bhns::set_x_rot ( )
inline

Sets the absolute coordinate X of the rotation axis [{ r_unit}].

Definition at line 215 of file bin_bhns.h.

References x_rot.

◆ set_y_rot()

double& Lorene::Bin_bhns::set_y_rot ( )
inline

Sets the absolute coordinate Y of the rotation axis [{ r_unit}].

Definition at line 218 of file bin_bhns.h.

References y_rot.

◆ shift_analytic()

◆ virial_bhns_surf()

double Lorene::Bin_bhns::virial_bhns_surf ( ) const

Estimates the relative error on the virial theorem $|1 - M_{ Komar} / M_{ ADM}|$.

Definition at line 2175 of file bin_bhns_global.C.

References mass_adm_bhns_surf(), mass_kom_bhns_surf(), and p_virial_bhns_surf.

◆ virial_bhns_vol()

double Lorene::Bin_bhns::virial_bhns_vol ( ) const

Estimates the relative error on the virial theorem $|1 - M_{ Komar} / M_{ ADM}|$.

Definition at line 2190 of file bin_bhns_global.C.

References p_virial_bhns_vol.

◆ xa_barycenter()

◆ ya_barycenter()

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  ost,
const Bin_bhns bibi 
)
friend

Display.

Definition at line 227 of file bin_bhns.C.

Member Data Documentation

◆ hole

Hole_bhns Lorene::Bin_bhns::hole
protected

Black hole.

Definition at line 72 of file bin_bhns.h.

◆ omega

double Lorene::Bin_bhns::omega
protected

Angular velocity with respect to an asymptotically inertial observer.

Definition at line 80 of file bin_bhns.h.

◆ p_angu_mom_bhns

Tbl* Lorene::Bin_bhns::p_angu_mom_bhns
mutableprotected

Total angular momentum of the system.

Definition at line 118 of file bin_bhns.h.

◆ p_line_mom_bhns

Tbl* Lorene::Bin_bhns::p_line_mom_bhns
mutableprotected

Total linear momentum of the system.

Definition at line 115 of file bin_bhns.h.

◆ p_mass_adm_bhns_surf

double* Lorene::Bin_bhns::p_mass_adm_bhns_surf
mutableprotected

Total ADM mass of the system calculated by the surface integral at infinity.

Definition at line 97 of file bin_bhns.h.

◆ p_mass_adm_bhns_vol

double* Lorene::Bin_bhns::p_mass_adm_bhns_vol
mutableprotected

Total ADM mass of the system calculated by the volume integral and the surface integral at the apparent horizon.

Definition at line 102 of file bin_bhns.h.

◆ p_mass_kom_bhns_surf

double* Lorene::Bin_bhns::p_mass_kom_bhns_surf
mutableprotected

Total Komar mass of the system calculated by the surface integral at infinity.

Definition at line 107 of file bin_bhns.h.

◆ p_mass_kom_bhns_vol

double* Lorene::Bin_bhns::p_mass_kom_bhns_vol
mutableprotected

Total Komar mass of the system calculated by the volume integral and the surface integral at the apparent horizon.

Definition at line 112 of file bin_bhns.h.

◆ p_omega_two_points

double* Lorene::Bin_bhns::p_omega_two_points
mutableprotected

Orbital angular velocity derived from another method.

Definition at line 137 of file bin_bhns.h.

◆ p_virial_bhns_surf

double* Lorene::Bin_bhns::p_virial_bhns_surf
mutableprotected

Virial theorem error calculated by the ADM mass and the Komar mass of the surface integral at infinity.

Definition at line 123 of file bin_bhns.h.

◆ p_virial_bhns_vol

double* Lorene::Bin_bhns::p_virial_bhns_vol
mutableprotected

Virial theorem error calculated by the ADM mass and the Komar mass of the volume integral.

Definition at line 128 of file bin_bhns.h.

◆ p_xa_barycenter

double* Lorene::Bin_bhns::p_xa_barycenter
mutableprotected

Absolute coordinate X of the barycenter of the baryon density.

Definition at line 131 of file bin_bhns.h.

◆ p_ya_barycenter

double* Lorene::Bin_bhns::p_ya_barycenter
mutableprotected

Absolute coordinate Y of the barycenter of the baryon density.

Definition at line 134 of file bin_bhns.h.

◆ ref_triad

const Base_vect_cart Lorene::Bin_bhns::ref_triad
protected

Cartesian triad of the absolute reference frame.

Definition at line 69 of file bin_bhns.h.

◆ separ

double Lorene::Bin_bhns::separ
protected

Absolute orbital separation between two centers of BH and NS.

Definition at line 83 of file bin_bhns.h.

◆ star

Star_bhns Lorene::Bin_bhns::star
protected

Neutron star.

Definition at line 75 of file bin_bhns.h.

◆ x_rot

double Lorene::Bin_bhns::x_rot
protected

Absolute X coordinate of the rotation axis.

Definition at line 86 of file bin_bhns.h.

◆ y_rot

double Lorene::Bin_bhns::y_rot
protected

Absolute Y coordinate of the rotation axis.

Definition at line 89 of file bin_bhns.h.


The documentation for this class was generated from the following files: