LORENE
Lorene::Bin_bhns_extr Class Reference

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. More...
 
 Bin_bhns_extr (const Bin_bhns_extr &)
 Copy constructor. More...
 
 Bin_bhns_extr (Map &mp, const Eos &eos, FILE *fich)
 Constructor from a file (see sauve(FILE*) ) More...
 
 ~Bin_bhns_extr ()
 Destructor. More...
 
void operator= (const Bin_bhns_extr &)
 Assignment to another Bin_bhns_extr. More...
 
Et_bin_bhns_extrset_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_mass_bh ()
 Sets the gravitational mass of BH [{ m_unit}]. More...
 
const Et_bin_bhns_extrget_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_mass_bh () const
 Returns the gravitational mass of BH [{ m_unit}]. More...
 
void sauve (FILE *) const
 Save in a file. More...
 
void display_poly (ostream &) const
 Display in polytropic units. More...
 
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. More...
 
double ya_barycenter_extr () const
 in the Kerr-Schild background metric More...
 
double mass_b_extr () const
 Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat. More...
 
void orbit_omega_ks (double fact_omeg_min, double fact_omeg_max)
 Computes the orbital angular velocity { omega} in the Kerr-Schild background metric. More...
 
void orbit_omega_cf (double fact_omeg_min, double fact_omeg_max)
 Computes the orbital angular velocity { omega} in the conformally flat background metric. More...
 
void analytical_omega ()
 Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian case) More...
 
void analytical_shift ()
 Sets some analytical template for the shift vector (via the members { w_shift} and { khi_shift} of { Etoile_bin}) More...
 

Private 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...
 
ostream & operator>> (ostream &) const
 Operator >> (function called by the operator <<) More...
 

Private Attributes

const Base_vect_cart ref_triad
 Cartesian triad of the absolute reference frame. More...
 
Et_bin_bhns_extr 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 mass_bh
 Gravitational mass of BH. More...
 
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. More...
 
double * p_ya_barycenter_extr
 Absolute coordinate Y of the barycenter of the baryon density in the Kerr-Schild background metric. More...
 
double * p_mass_b_extr
 Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat one. More...
 

Friends

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

Detailed Description

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.

Constructor & Destructor Documentation

◆ Bin_bhns_extr() [1/3]

Lorene::Bin_bhns_extr::Bin_bhns_extr ( Map mp,
int  nzet,
const Eos eos,
bool  irrot,
bool  relat,
bool  kerrs,
bool  multi 
)

Standard constructor.

Parameters
mpMapping on which the neutron star will be defined
nzetNumber of domains occupied by the neutron star
eosEquation of state of the neutron star
irrotshould be { true} if NS is irrotational, { false} if NS is corotating
relatshould be { true} for a relativistic configuration, { false} for a Newtonian one
kerrsshould be { true} for the Kerr-Schild background metric, { false} for the conformally flat one
multishould be { true} for the multipole falloff boundary condition, { true} for the $1/r$ one

Definition at line 88 of file bin_bhns_extr.C.

References mass_bh, omega, separ, and set_der_0x0().

◆ Bin_bhns_extr() [2/3]

Lorene::Bin_bhns_extr::Bin_bhns_extr ( const Bin_bhns_extr bibi)

Copy constructor.

Definition at line 105 of file bin_bhns_extr.C.

References set_der_0x0().

◆ Bin_bhns_extr() [3/3]

Lorene::Bin_bhns_extr::Bin_bhns_extr ( Map mp,
const Eos eos,
FILE *  fich 
)

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

Parameters
mpMapping on which the neutron star will be defined
eosEquation of state of the neutron star

Definition at line 120 of file bin_bhns_extr.C.

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

◆ ~Bin_bhns_extr()

Lorene::Bin_bhns_extr::~Bin_bhns_extr ( )

Destructor.

Definition at line 139 of file bin_bhns_extr.C.

References del_deriv().

Member Function Documentation

◆ analytical_omega()

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 58 of file bin_bhns_extr_omegaana.C.

References del_deriv(), Lorene::Etoile_bin::is_irrotational(), Lorene::Etoile::is_relativistic(), mass_bh, omega, Lorene::pow(), Lorene::Etoile::ray_eq(), separ, Lorene::sqrt(), and star.

◆ analytical_shift()

◆ del_deriv()

void Lorene::Bin_bhns_extr::del_deriv ( ) const
private

Deletes all the derived quantities.

Definition at line 148 of file bin_bhns_extr.C.

References p_mass_b_extr, p_xa_barycenter_extr, p_ya_barycenter_extr, and set_der_0x0().

◆ display_poly()

◆ get_mass_bh()

double Lorene::Bin_bhns_extr::get_mass_bh ( ) const
inline

Returns the gravitational mass of BH [{ m_unit}].

Definition at line 178 of file bin_bhns_extr.h.

References mass_bh.

◆ get_ns()

const Et_bin_bhns_extr& Lorene::Bin_bhns_extr::get_ns ( ) const
inline

Returns a reference to the neutron star.

Definition at line 166 of file bin_bhns_extr.h.

References star.

◆ get_omega()

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

Returns the orbital angular velocity [{ f_unit}].

Definition at line 170 of file bin_bhns_extr.h.

References omega.

◆ get_separ()

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

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

Definition at line 175 of file bin_bhns_extr.h.

References separ.

◆ mass_b_extr()

◆ operator=()

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

Assignment to another Bin_bhns_extr.

Definition at line 172 of file bin_bhns_extr.C.

References del_deriv(), mass_bh, omega, ref_triad, separ, and star.

◆ operator>>()

◆ orbit_omega_cf()

void Lorene::Bin_bhns_extr::orbit_omega_cf ( double  fact_omeg_min,
double  fact_omeg_max 
)

Computes the orbital angular velocity { omega} in the conformally flat background metric.

Parameters
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})
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 325 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::Etoile::get_a_car(), Lorene::Map::get_bvect_cart(), Lorene::Etoile::get_d_logn_auto_div(), Lorene::Etoile_bin::get_loggam(), Lorene::Etoile::get_logn_auto_regu(), Lorene::Etoile_bin::get_logn_comp(), Lorene::Etoile::get_mp(), Lorene::Etoile::get_nnn(), Lorene::Map::get_rot_phi(), Lorene::Etoile::get_shift(), Lorene::Tbl::get_taille(), Lorene::Tenseur::get_triad(), Lorene::Etoile::is_relativistic(), omega, ref_triad, star, Lorene::zero_list(), and Lorene::zerosec_b().

◆ orbit_omega_ks()

void Lorene::Bin_bhns_extr::orbit_omega_ks ( double  fact_omeg_min,
double  fact_omeg_max 
)

Computes the orbital angular velocity { omega} in the Kerr-Schild background metric.

Parameters
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})
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 73 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::Etoile::get_a_car(), Lorene::Map::get_bvect_cart(), Lorene::Etoile::get_d_logn_auto_div(), Lorene::Etoile_bin::get_d_logn_comp(), Lorene::Etoile_bin::get_loggam(), Lorene::Etoile::get_logn_auto_regu(), Lorene::Etoile::get_mp(), Lorene::Etoile::get_nnn(), Lorene::Map::get_rot_phi(), Lorene::Etoile::get_shift(), Lorene::Tbl::get_taille(), Lorene::Tenseur::get_triad(), Lorene::Etoile::is_relativistic(), mass_bh, omega, ref_triad, separ, star, Lorene::zero_list(), and Lorene::zerosec_b().

◆ sauve()

void Lorene::Bin_bhns_extr::sauve ( FILE *  fich) const

Save in a file.

Definition at line 194 of file bin_bhns_extr.C.

References Lorene::fwrite_be(), mass_bh, omega, Lorene::Et_bin_bhns_extr::sauve(), separ, and star.

◆ set_der_0x0()

void Lorene::Bin_bhns_extr::set_der_0x0 ( ) const
private

Sets to 0x0 all the pointers on derived quantities.

Definition at line 158 of file bin_bhns_extr.C.

References p_mass_b_extr, p_xa_barycenter_extr, and p_ya_barycenter_extr.

◆ set_mass_bh()

double& Lorene::Bin_bhns_extr::set_mass_bh ( )
inline

Sets the gravitational mass of BH [{ m_unit}].

Definition at line 160 of file bin_bhns_extr.h.

References mass_bh.

◆ set_ns()

Et_bin_bhns_extr& Lorene::Bin_bhns_extr::set_ns ( )
inline

Read/write of the neutron star.

Definition at line 149 of file bin_bhns_extr.h.

References del_deriv(), and star.

◆ set_omega()

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

Sets the orbital angular velocity [{ f_unit}].

Definition at line 154 of file bin_bhns_extr.h.

References omega.

◆ set_separ()

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

Sets the orbital separation [{ r_unit}].

Definition at line 157 of file bin_bhns_extr.h.

References separ.

◆ xa_barycenter_extr()

◆ ya_barycenter_extr()

Friends And Related Function Documentation

◆ operator<<

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

Display.

Definition at line 206 of file bin_bhns_extr.C.

Member Data Documentation

◆ mass_bh

double Lorene::Bin_bhns_extr::mass_bh
private

Gravitational mass of BH.

Definition at line 77 of file bin_bhns_extr.h.

◆ omega

double Lorene::Bin_bhns_extr::omega
private

Angular velocity with respect to an asymptotically inertial observer.

Definition at line 71 of file bin_bhns_extr.h.

◆ p_mass_b_extr

double* Lorene::Bin_bhns_extr::p_mass_b_extr
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.

◆ p_xa_barycenter_extr

double* Lorene::Bin_bhns_extr::p_xa_barycenter_extr
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.

◆ p_ya_barycenter_extr

double* Lorene::Bin_bhns_extr::p_ya_barycenter_extr
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.

◆ ref_triad

const Base_vect_cart Lorene::Bin_bhns_extr::ref_triad
private

Cartesian triad of the absolute reference frame.

Definition at line 63 of file bin_bhns_extr.h.

◆ separ

double Lorene::Bin_bhns_extr::separ
private

Absolute orbital separation between two centers of BH and NS.

Definition at line 74 of file bin_bhns_extr.h.

◆ star

Et_bin_bhns_extr Lorene::Bin_bhns_extr::star
private

Neutron star.

Definition at line 66 of file bin_bhns_extr.h.


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