LORENE
|
Binary systems in eXtended Conformal Thin Sandwich formulation. More...
#include <binary_xcts.h>
Public Member Functions | |
Binary_xcts (Map &mp1, int nzet1, const Eos &eos1, int irrot1, Map &mp2, int nzet2, const Eos &eos2, int irrot2) | |
Standard constructor. More... | |
Binary_xcts (const Binary_xcts &) | |
Copy constructor. More... | |
Binary_xcts (Map &mp1, const Eos &eos1, Map &mp2, const Eos &eos2, FILE *fich) | |
Constructor from a file (see sauve(FILE* ) ). More... | |
~Binary_xcts () | |
Destructor. More... | |
void | operator= (const Binary_xcts &) |
Assignment to another Binary_xcts . More... | |
Star_bin_xcts & | set (int i) |
Read/write of the star no. i. More... | |
double & | set_omega () |
Sets the orbital angular velocity [f_unit ]. More... | |
double & | set_x_axe () |
Sets the absolute coordinate X of the rotation axis [r_unit ]. More... | |
const Star_bin_xcts & | operator() (int i) const |
Returns a reference to the star no. i. More... | |
double | get_omega () const |
Returns the orbital angular velocity [f_unit ]. More... | |
double | get_x_axe () const |
Returns the absolute coordinate X of the rotation axis [r_unit ]. More... | |
double | separation () const |
Returns the coordinate separation of the two stellar centers [r_unit ]. More... | |
void | sauve (FILE *) const |
Save in a file. More... | |
void | display_poly (ostream &) const |
Display in polytropic units. More... | |
void | write_global (ostream &) const |
Write global quantities in a formatted file. More... | |
double | mass_adm () const |
Total ADM mass. More... | |
double | mass_adm_vol () const |
Total ADM mass (computed by a volume integral) More... | |
double | mass_kom () const |
Total Komar mass. More... | |
double | mass_kom_vol () const |
Total Komar mass (computed by a volume integral) More... | |
const Tbl & | angu_mom () const |
Total angular momentum. More... | |
const Tbl & | lin_mom () const |
Total linear momentum. More... | |
double | total_ener () const |
Total energy (excluding the rest mass energy). More... | |
double | virial () const |
Estimates the relative error on the virial theorem. More... | |
double | virial_vol () const |
Estimates the relative error on the virial theorem (volume version) More... | |
double | ham_constr () const |
Estimates the relative error on the Hamiltonian constraint. More... | |
const Tbl & | mom_constr () const |
Estimates the relative error on the momentum constraint. More... | |
void | orbit (double fact_omeg_min, double fact_omeg_max, double &xgg1, double &xgg2) |
Computes the orbital angular velocity omega and the position of the rotation axis x_axe . 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 the two Star_bin . 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 | |
Star_bin_xcts | star1 |
First star of the system. More... | |
Star_bin_xcts | star2 |
Second star of the system. More... | |
Star_bin_xcts * | et [2] |
Array of the two stars (to perform loops on the stars): et [0] contains the address of star1 and et [1] that of star2 . More... | |
double | omega |
Angular velocity with respect to an asymptotically inertial observer. More... | |
double | x_axe |
Absolute X coordinate of the rotation axis. More... | |
double * | p_mass_adm |
Total ADM mass of the system. More... | |
double * | p_mass_kom |
Total Komar mass of the system. More... | |
Tbl * | p_angu_mom |
Total angular momentum of the system. More... | |
Tbl * | p_lin_mom |
Total linear momentum of the system. More... | |
double * | p_total_ener |
Total energy of the system. More... | |
double * | p_virial |
Virial theorem error. More... | |
double * | p_virial_vol |
Virial theorem error (volume version) More... | |
Private Member Functions | |
ostream & | operator>> (ostream &) const |
Operator >> (function called by the operator <<). More... | |
Friends | |
ostream & | operator<< (ostream &, const Binary_xcts &) |
Display. More... | |
Binary systems in eXtended Conformal Thin Sandwich formulation.
*** UNDER DEVELOPMENT ***()
Definition at line 59 of file binary_xcts.h.
Lorene::Binary_xcts::Binary_xcts | ( | Map & | mp1, |
int | nzet1, | ||
const Eos & | eos1, | ||
int | irrot1, | ||
Map & | mp2, | ||
int | nzet2, | ||
const Eos & | eos2, | ||
int | irrot2 | ||
) |
Standard constructor.
mp1 | Mapping on which star1 will be defined |
nzet1 | Number of domains occupied by star1 |
eos1 | Equation of state of star1 |
irrot1 | should be true if star1 is irrotational, false if star1 is corotating |
mp2 | Mapping on which star2 will be defined |
nzet2 | Number of domains occupied by star2 |
eos2 | Equation of state of star2 |
irrot2 | should be true if star2 is irrotational, false if star2 is corotating |
Definition at line 75 of file binary_xcts.C.
References et, omega, set_der_0x0(), star1, star2, and x_axe.
Lorene::Binary_xcts::Binary_xcts | ( | const Binary_xcts & | bibi | ) |
Copy constructor.
Definition at line 99 of file binary_xcts.C.
References et, set_der_0x0(), star1, and star2.
Lorene::Binary_xcts::Binary_xcts | ( | Map & | mp1, |
const Eos & | eos1, | ||
Map & | mp2, | ||
const Eos & | eos2, | ||
FILE * | fich | ||
) |
Constructor from a file (see sauve(FILE* )
).
mp1 | Mapping on which star1 will be defined |
eos1 | Equation of state of star1 |
mp2 | Mapping on which star2 will be defined |
eos2 | Equation of state of star2 |
fich | input file (must have been created by the function sauve ) |
Definition at line 114 of file binary_xcts.C.
References et, Lorene::fread_be(), omega, set_der_0x0(), star1, star2, and x_axe.
Lorene::Binary_xcts::~Binary_xcts | ( | ) |
void Lorene::Binary_xcts::analytical_omega | ( | ) |
Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian case)
Definition at line 53 of file binary_omega_ana_xcts.C.
References del_deriv(), Lorene::Star_bin_xcts::is_irrotational(), Lorene::Star_bin_xcts::mass_g(), omega, Lorene::pow(), Lorene::Star::ray_eq(), separation(), Lorene::sqrt(), star1, and star2.
void Lorene::Binary_xcts::analytical_shift | ( | ) |
Sets some analytical template for the shift vector (via the members w_shift
and khi_shift
of the two Star_bin
.
Definition at line 56 of file binary_anashift_xcts.C.
References Lorene::Scalar::annule(), Lorene::Tenseur::dec_dzpuis(), Lorene::Tensor::dec_dzpuis(), Lorene::Tensor::derive_con(), et, Lorene::Map::get_bvect_cart(), Lorene::Star_bin_xcts::get_flat(), Lorene::Map::get_mg(), Lorene::Star::get_mp(), Lorene::Star::get_nzet(), Lorene::Mg3d::get_nzone(), Lorene::Star_bin_xcts::mass_g(), omega, Lorene::Map::r, Lorene::Star::ray_eq(), separation(), Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Tenseur::set(), Lorene::Star_bin_xcts::set_beta_auto(), Lorene::Tenseur::set_etat_qcq(), Lorene::Map::sinp, Lorene::Map::sint, Lorene::skxk(), Lorene::Vector::std_spectral_base(), Lorene::Scalar::std_spectral_base(), and Lorene::Map::y.
const Tbl & Lorene::Binary_xcts::angu_mom | ( | ) | const |
Total angular momentum.
Tbl
of size 3, according to angu_mom()
(0) = , angu_mom()
(1) = , angu_mom()
(2) = . Definition at line 216 of file binary_global_xcts.C.
References Lorene::Tbl::annule_hard(), Lorene::Vector::change_triad(), Lorene::contract(), et, Lorene::Map::get_bvect_cart(), Lorene::Star::get_ener_euler(), Lorene::Star_bin_xcts::get_flat(), Lorene::Star::get_mp(), Lorene::Star::get_press(), Lorene::Scalar::integrale(), p_angu_mom, Lorene::pow(), Lorene::Tbl::set(), Lorene::Vector::set(), Lorene::Vector::std_spectral_base(), Lorene::Tensor::up_down(), Lorene::Map::xa, and Lorene::Map::ya.
|
protected |
Deletes all the derived quantities.
Definition at line 149 of file binary_xcts.C.
References p_angu_mom, p_lin_mom, p_mass_adm, p_mass_kom, p_total_ener, p_virial, p_virial_vol, and set_der_0x0().
void Lorene::Binary_xcts::display_poly | ( | ostream & | ost | ) | const |
Display in polytropic units.
Definition at line 249 of file binary_xcts.C.
References angu_mom(), Lorene::Star::get_eos(), Lorene::Eos_poly::get_gam(), Lorene::Eos_poly::get_kap(), mass_adm(), Lorene::Star_bin_xcts::mass_b(), mass_kom(), omega, Lorene::pow(), Lorene::Star::ray_eq(), Lorene::Star::ray_eq_pi(), separation(), Lorene::sqrt(), star1, star2, and Lorene::Star_bin_xcts::xa_barycenter().
|
inline |
Returns the orbital angular velocity [f_unit
].
Definition at line 190 of file binary_xcts.h.
References omega.
|
inline |
Returns the absolute coordinate X of the rotation axis [r_unit
].
Definition at line 193 of file binary_xcts.h.
References x_axe.
double Lorene::Binary_xcts::ham_constr | ( | ) | const |
Estimates the relative error on the Hamiltonian constraint.
const Tbl & Lorene::Binary_xcts::lin_mom | ( | ) | const |
Total linear momentum.
Definition at line 271 of file binary_global_xcts.C.
References Lorene::Tbl::annule_hard(), Lorene::Vector::change_triad(), et, Lorene::Star::get_ener_euler(), Lorene::Star::get_press(), p_lin_mom, Lorene::pow(), Lorene::Tbl::set(), and Lorene::Vector::std_spectral_base().
double Lorene::Binary_xcts::mass_adm | ( | ) | const |
Total ADM mass.
Definition at line 85 of file binary_global_xcts.C.
References Lorene::Vector::change_triad(), et, Lorene::Map::get_bvect_spher(), Lorene::Star_bin_xcts::get_flat(), Lorene::Map_af::integrale_surface_infini(), and p_mass_adm.
double Lorene::Binary_xcts::mass_adm_vol | ( | ) | const |
Total ADM mass (computed by a volume integral)
Definition at line 116 of file binary_global_xcts.C.
References et, Lorene::Star::get_ener_euler(), Lorene::Star_bin_xcts::get_hacar_auto(), Lorene::Star_bin_xcts::get_hacar_comp(), Lorene::Scalar::integrale(), Lorene::pow(), and Lorene::Scalar::std_spectral_base().
double Lorene::Binary_xcts::mass_kom | ( | ) | const |
Total Komar mass.
Definition at line 152 of file binary_global_xcts.C.
References Lorene::Vector::change_triad(), Lorene::Scalar::derive_con(), et, Lorene::Map::get_bvect_spher(), Lorene::Star_bin_xcts::get_flat(), Lorene::Star::get_logn(), Lorene::Map_af::integrale_surface_infini(), and p_mass_kom.
double Lorene::Binary_xcts::mass_kom_vol | ( | ) | const |
Total Komar mass (computed by a volume integral)
Definition at line 179 of file binary_global_xcts.C.
References et, Lorene::Star_bin_xcts::get_chi(), Lorene::Star::get_ener_euler(), Lorene::Star_bin_xcts::get_hacar_auto(), Lorene::Star_bin_xcts::get_hacar_comp(), Lorene::Star_bin_xcts::get_Psi(), Lorene::Star_bin_xcts::get_psi4(), Lorene::Star::get_s_euler(), Lorene::Scalar::integrale(), Lorene::pow(), and Lorene::Scalar::std_spectral_base().
const Tbl& Lorene::Binary_xcts::mom_constr | ( | ) | const |
Estimates the relative error on the momentum constraint.
|
inline |
void Lorene::Binary_xcts::operator= | ( | const Binary_xcts & | bibi | ) |
Assignment to another Binary_xcts
.
Definition at line 185 of file binary_xcts.C.
References del_deriv(), omega, star1, star2, and x_axe.
|
private |
Operator >> (function called by the operator <<).
Definition at line 222 of file binary_xcts.C.
References omega, separation(), star1, star2, and x_axe.
void Lorene::Binary_xcts::orbit | ( | double | fact_omeg_min, |
double | fact_omeg_max, | ||
double & | xgg1, | ||
double & | xgg2 | ||
) |
Computes the orbital angular velocity omega
and the position of the rotation axis x_axe
.
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) |
xgg1 | [output] : x coordinate (relative to star 1 mapping) of the `‘center of mass’' of star 1 |
xgg2 | [output] : x coordinate (relative to star 2 mapping) of the `‘center of mass’' of star 2 |
Definition at line 100 of file binary_orbit_xcts.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Vector::change_triad(), Lorene::contract(), Lorene::Scalar::dsdx(), et, Lorene::Star::get_beta(), Lorene::Map::get_bvect_cart(), Lorene::Star_bin_xcts::get_chi_auto(), Lorene::Star_bin_xcts::get_chi_comp(), Lorene::Star::get_eos(), Lorene::Star_bin_xcts::get_flat(), Lorene::Star_bin_xcts::get_loggam(), Lorene::Star::get_mp(), Lorene::Star_bin_xcts::get_Psi_auto(), Lorene::Star_bin_xcts::get_Psi_comp(), Lorene::Map::get_rot_phi(), Lorene::Tbl::get_taille(), Lorene::log(), omega, Lorene::pow(), save_profile(), Lorene::Scalar::std_spectral_base(), Lorene::Tensor::up_down(), Lorene::Scalar::val_grid_point(), x_axe, Lorene::Star_bin_xcts::xa_barycenter(), Lorene::zero_list(), Lorene::zerosec(), and Lorene::zerosec_b().
void Lorene::Binary_xcts::sauve | ( | FILE * | fich | ) | const |
Save in a file.
Definition at line 203 of file binary_xcts.C.
References Lorene::fwrite_be(), omega, Lorene::Star_bin_xcts::sauve(), star1, star2, and x_axe.
double Lorene::Binary_xcts::separation | ( | ) | const |
Returns the coordinate separation of the two stellar centers [r_unit
].
Definition at line 437 of file binary_xcts.C.
References Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), Lorene::Map::get_ori_z(), Lorene::Star::mp, Lorene::sqrt(), star1, and star2.
|
inline |
Read/write of the star no. i.
Definition at line 167 of file binary_xcts.h.
References del_deriv(), and et.
|
protected |
Sets to 0x0
all the pointers on derived quantities.
Definition at line 165 of file binary_xcts.C.
References p_angu_mom, p_lin_mom, p_mass_adm, p_mass_kom, p_total_ener, p_virial, and p_virial_vol.
|
inline |
Sets the orbital angular velocity [f_unit
].
Definition at line 175 of file binary_xcts.h.
References omega.
|
inline |
Sets the absolute coordinate X of the rotation axis [r_unit
].
Definition at line 178 of file binary_xcts.h.
References x_axe.
double Lorene::Binary_xcts::total_ener | ( | ) | const |
Total energy (excluding the rest mass energy).
In the Newtonian case, it is defined as the sum of kinetic, internal, and gravitational potential energies.
In the relativistic case, it is defined as .
Definition at line 310 of file binary_global_xcts.C.
References mass_adm(), Lorene::Star_bin_xcts::mass_b(), p_total_ener, star1, and star2.
double Lorene::Binary_xcts::virial | ( | ) | const |
Estimates the relative error on the virial theorem.
Definition at line 329 of file binary_global_xcts.C.
References mass_adm(), mass_kom(), and p_virial.
double Lorene::Binary_xcts::virial_vol | ( | ) | const |
Estimates the relative error on the virial theorem (volume version)
Definition at line 348 of file binary_global_xcts.C.
References mass_adm_vol(), mass_kom_vol(), and p_virial_vol.
void Lorene::Binary_xcts::write_global | ( | ostream & | ost | ) | const |
Write global quantities in a formatted file.
This file can be read by an external program.
Definition at line 299 of file binary_xcts.C.
References angu_mom(), Lorene::Star::ent, Lorene::Star::get_ener(), Lorene::Star::get_ent(), Lorene::Star::get_eos(), Lorene::Eos_poly::get_gam(), Lorene::Eos_poly::get_kap(), Lorene::Map::get_mg(), Lorene::Star::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), mass_adm(), mass_adm_vol(), Lorene::Star_bin_xcts::mass_b(), mass_kom(), mass_kom_vol(), omega, Lorene::pow(), Lorene::Star::ray_eq(), Lorene::Star::ray_eq_pis2(), Lorene::Star::ray_pole(), separation(), Lorene::sqrt(), star1, star2, Lorene::Scalar::val_grid_point(), Lorene::Map::val_r(), virial(), virial_vol(), and Lorene::Star_bin_xcts::xa_barycenter().
|
friend |
Display.
Definition at line 215 of file binary_xcts.C.
|
protected |
Array of the two stars (to perform loops on the stars): et
[0] contains the address of star1
and et
[1] that of star2
.
Definition at line 75 of file binary_xcts.h.
|
protected |
Angular velocity with respect to an asymptotically inertial observer.
Definition at line 80 of file binary_xcts.h.
|
mutableprotected |
Total angular momentum of the system.
Definition at line 97 of file binary_xcts.h.
|
mutableprotected |
Total linear momentum of the system.
Definition at line 100 of file binary_xcts.h.
|
mutableprotected |
Total ADM mass of the system.
Definition at line 91 of file binary_xcts.h.
|
mutableprotected |
Total Komar mass of the system.
Definition at line 94 of file binary_xcts.h.
|
mutableprotected |
Total energy of the system.
Definition at line 103 of file binary_xcts.h.
|
mutableprotected |
Virial theorem error.
Definition at line 106 of file binary_xcts.h.
|
mutableprotected |
Virial theorem error (volume version)
Definition at line 109 of file binary_xcts.h.
|
protected |
First star of the system.
Definition at line 66 of file binary_xcts.h.
|
protected |
Second star of the system.
Definition at line 69 of file binary_xcts.h.
|
protected |
Absolute X coordinate of the rotation axis.
Definition at line 84 of file binary_xcts.h.