Binaire Class Reference

Binary systems. More...

#include <binaire.h>

List of all members.

Public Member Functions

 Binaire (Map &mp1, int nzet1, const Eos &eos1, int irrot1, Map &mp2, int nzet2, const Eos &eos2, int irrot2, int relat)
 Standard constructor.
 Binaire (const Binaire &)
 Binaire (Map &mp1, const Eos &eos1, Map &mp2, const Eos &eos2, FILE *fich)
 Copy constructor.
void operator= (const Binaire &)
 Assignment to another { Binaire}.
Etoile_binset (int i)
 Read/write of the star no. i.
double & set_omega ()
 Sets the orbital angular velocity [{ f}].
double & set_x_axe ()
 Sets the absolute coordinate X of the rotation axis [{ r}].
const Etoile_binoperator() (int i) const
 Returns a reference to the star no. i.
double get_omega () const
 Returns the orbital angular velocity [{ f}].
double get_x_axe () const
 Returns the absolute coordinate X of the rotation axis [{ r}].
double separation () const
 Returns the coordinate separation of the two stellar centers [{ r}].
void sauve (FILE *) const
void display_poly (ostream &) const
 Display in polytropic units.
void write_global (ostream &) const
 Write global quantities in a formatted file.
double mass_adm () const
 Total ADM mass.
double mass_kom () const
 Total Komar mass.
const Tblangu_mom () const
 Total angular momentum.
double total_ener () const
 Total energy (excluding the rest mass energy).
double virial () const
 Estimates the relative error on the virial theorem (for a relativistic one, it returns $|1 - M_{ Komar} / M_{ ADM}|$).
double virial_gb () const
 Estimates the relative error on the virial theorem calculated by E.Gourgoulhon and S.Bonazzola (Class.
double virial_fus () const
 Estimates the relative error on the virial theorem calculated by J.L.Friedman, K.Uryu, and M.Shibata (PRD accepted, gr-qc/0108070).
double ham_constr () const
 Estimates the relative error on the Hamiltonian constraint equation by comparing $ A$ with {equation} -4 A^2 E - {A^2 4} K_{ij} K^{ij} - {1 2} {}_i A {}^i A {equation}.
const Tblmom_constr () const
 Estimates the relative error on the momentum constraint equation by comparing ${}_j K^{ij}$ with {equation} 8 J^i - 5 K^{ij} {}_j A {equation}.
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}.
void orbit_eqmass (double fact_omeg_min, double fact_omeg_max, double mass1, double mass2, double &xgg1, double &xgg2)
 Computes the orbital angular velocity { omega} and the position of the rotation axis { x}.
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 { w} and { khi} of the two { Etoile}.

Private Member Functions

void del_deriv () const
 Destructor.
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.
Etoile_bin star1
 First star of the system.
Etoile_bin star2
 Second star of the system.
Etoile_binet [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}.
double omega
 Angular velocity with respect to an asymptotically inertial observer.
double x_axe
 Absolute X coordinate of the rotation axis.
double * p_mass_adm
 Total ADM mass of the system.
double * p_mass_kom
 Total Komar mass of the system.
Tblp_angu_mom
 Total angular momentum of the system.
double * p_total_ener
 Total energy of the system.
double * p_virial
 Virial theorem error.
double * p_virial_gb
 Virial theorem error by E.Gourgoulhon and S.Bonazzola.
double * p_virial_fus
 Virial theorem error by J.L.Friedman, K.Uryu, and M.Shibata.
double * p_ham_constr
 Relative error on the Hamiltonian constraint.
Tblp_mom_constr
 Relative error on the momentum constraint.

Friends

ostream & operator<< (ostream &, const Binaire &)
 Save in a file.

Detailed Description

Binary systems.

Version:
$Id: binaire.h,v 1.5 2009/06/18 18:42:13 k_taniguchi Exp $#

Definition at line 100 of file binaire.h.


Constructor & Destructor Documentation

Binaire::Binaire ( Map mp1,
int  nzet1,
const Eos eos1,
int  irrot1,
Map mp2,
int  nzet2,
const Eos eos2,
int  irrot2,
int  relat 
)

Standard constructor.

Parameters:
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
relat should be { true} for a relativistic configuration, { false} for a Newtonian one

Definition at line 105 of file binaire.C.

References et, omega, set_der_0x0(), star1, star2, and x_axe.

Binaire::Binaire ( Map mp1,
const Eos eos1,
Map mp2,
const Eos eos2,
FILE *  fich 
)

Copy constructor.

Constructor from a file (see { sauve(FILE* )}).

Parameters:
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 141 of file binaire.C.

References et, fread_be(), omega, set_der_0x0(), star1, star2, and x_axe.


Member Function Documentation

void Binaire::analytical_omega (  ) 

Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian case).

Definition at line 59 of file binaire_omega_ana.C.

References del_deriv(), Etoile_bin::is_irrotational(), Etoile::is_relativistic(), Etoile_bin::mass_g(), omega, pow(), Etoile::ray_eq(), separation(), sqrt(), star1, and star2.

void Binaire::analytical_shift (  ) 
const Tbl & Binaire::angu_mom (  )  const

Total angular momentum.

Returns:
1-D { Tbl} of size 3, according to \ { angu()(0)} = $J^X$, \ { angu()(1)} = $J^Y$, \ { angu()(2)} = $J^Z$.

Definition at line 194 of file binaire_global.C.

References Tbl::annule_hard(), et, Etoile::get_a_car(), Etoile::get_ener_euler(), Etoile::get_mp(), Etoile::get_nbar(), Etoile::get_press(), Etoile::get_u_euler(), Cmp::integrale(), p_angu_mom, pow(), Tbl::set(), Cmp::std_base_scal(), Cmp::va, Map::xa, Map::ya, and Map::za.

void Binaire::del_deriv (  )  const [private]

Destructor.

Deletes all the derived quantities

Definition at line 173 of file binaire.C.

References p_angu_mom, p_ham_constr, p_mass_adm, p_mass_kom, p_mom_constr, p_total_ener, p_virial, p_virial_fus, p_virial_gb, and set_der_0x0().

void Binaire::display_poly ( ostream &  ost  )  const
double Binaire::get_omega (  )  const [inline]

Returns the orbital angular velocity [{ f}].

Definition at line 240 of file binaire.h.

References omega.

double Binaire::get_x_axe (  )  const [inline]

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

Definition at line 243 of file binaire.h.

References x_axe.

double Binaire::ham_constr (  )  const

Estimates the relative error on the Hamiltonian constraint equation by comparing $ A$ with {equation} -4 A^2 E - {A^2 4} K_{ij} K^{ij} - {1 2} {}_i A {}^i A {equation}.

Definition at line 65 of file binaire_constr.C.

References abs(), diffrel(), et, flat_scalar_prod(), Etoile::get_a_car(), Etoile_bin::get_akcar_auto(), Etoile_bin::get_akcar_comp(), Etoile::get_beta_auto(), Etoile_bin::get_d_beta_auto(), Etoile_bin::get_d_beta_comp(), Etoile_bin::get_d_logn_auto(), Etoile_bin::get_d_logn_comp(), Etoile::get_ener_euler(), Etoile::get_logn_auto(), Etoile::get_mp(), max(), p_ham_constr, star1, and star2.

double Binaire::mass_adm (  )  const
double Binaire::mass_kom (  )  const
const Tbl & Binaire::mom_constr (  )  const
const Etoile_bin& Binaire::operator() ( int  i  )  const [inline]

Returns a reference to the star no. i.

Definition at line 235 of file binaire.h.

References et.

void Binaire::operator= ( const Binaire bibi  ) 

Assignment to another { Binaire}.

Definition at line 213 of file binaire.C.

References del_deriv(), omega, ref_triad, star1, star2, and x_axe.

ostream & Binaire::operator>> ( ostream &  ost  )  const [private]

Operator >> (function called by the operator <<).

Definition at line 253 of file binaire.C.

References omega, separation(), star1, star2, and x_axe.

void Binaire::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}.

Parameters:
fact_omeg_min [input] : determines the lower bound of the interval { [omega, omega]} in which { omega} is searched by { omega = fact * omega}, where { omega} is the previous value of the angular velocity (typical value : { fact = 0.5})
fact_omeg_max [input] : determines the higher bound of the interval { [omega, omega]} in which { omega} is searched by { omega = fact * omega}, where { omega} is the previous value of the angular velocity. (typical value : { fact = 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 115 of file binaire_orbite.C.

References Param::add_double(), Param::add_int(), Tenseur::change_triad(), Cmp::dsdx(), et, flat_scalar_prod(), Etoile::get_a_car(), Map::get_bvect_cart(), Etoile::get_d_logn_auto_div(), Etoile::get_ent(), Etoile::get_eos(), Etoile_bin::get_loggam(), Etoile::get_logn_auto_regu(), Etoile_bin::get_logn_comp(), Etoile::get_mp(), Etoile::get_nnn(), Map::get_rot_phi(), Etoile::get_shift(), Tbl::get_taille(), Tenseur::get_triad(), Etoile::is_relativistic(), omega, ref_triad, save_profile(), x_axe, Etoile_bin::xa_barycenter(), zero_list(), zerosec(), and zerosec_b().

void Binaire::orbit_eqmass ( double  fact_omeg_min,
double  fact_omeg_max,
double  mass1,
double  mass2,
double &  xgg1,
double &  xgg2 
)

Computes the orbital angular velocity { omega} and the position of the rotation axis { x}.

Parameters:
fact_omeg_min [input] : determines the lower bound of the interval { [omega, omega]} in which { omega} is searched by { omega = fact * omega}, where { omega} is the previous value of the angular velocity (typical value : { fact = 0.5})
fact_omeg_max [input] : determines the higher bound of the interval { [omega, omega]} in which { omega} is searched by { omega = fact * omega}, where { omega} is the previous value of the angular velocity. (typical value : { fact = 1.5})
mass1 [input] : baryon rest mass of NS1
mass2 [input] : baryon rest mass of NS2
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 457 of file binaire_orbite.C.

References Param::add_double(), Param::add_int(), Tenseur::change_triad(), Cmp::dsdx(), et, flat_scalar_prod(), Etoile::get_a_car(), Map::get_bvect_cart(), Etoile::get_d_logn_auto_div(), Etoile::get_eos(), Etoile_bin::get_loggam(), Etoile::get_logn_auto_regu(), Etoile_bin::get_logn_comp(), Etoile::get_mp(), Etoile::get_nnn(), Map::get_rot_phi(), Etoile::get_shift(), Tbl::get_taille(), Tenseur::get_triad(), Etoile::is_relativistic(), omega, ref_triad, x_axe, Etoile_bin::xa_barycenter(), zero_list(), zerosec(), and zerosec_b().

double Binaire::separation (  )  const

Returns the coordinate separation of the two stellar centers [{ r}].

Definition at line 453 of file binaire.C.

References Etoile::get_mp(), Map::get_ori_x(), Map::get_ori_y(), Map::get_ori_z(), sqrt(), star1, and star2.

Etoile_bin& Binaire::set ( int  i  )  [inline]

Read/write of the star no. i.

Definition at line 219 of file binaire.h.

References del_deriv(), and et.

void Binaire::set_der_0x0 (  )  const [private]

Sets to { 0x0} all the pointers on derived quantities.

Definition at line 191 of file binaire.C.

References p_angu_mom, p_ham_constr, p_mass_adm, p_mass_kom, p_mom_constr, p_total_ener, p_virial, p_virial_fus, and p_virial_gb.

double& Binaire::set_omega (  )  [inline]

Sets the orbital angular velocity [{ f}].

Definition at line 225 of file binaire.h.

References omega.

double& Binaire::set_x_axe (  )  [inline]

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

Definition at line 228 of file binaire.h.

References x_axe.

double Binaire::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 $M_{ ADM} - M_{ bar,1} - M_{ bar,2}$.

Definition at line 285 of file binaire_global.C.

References et, flat_scalar_prod(), Etoile::get_ener(), Etoile::get_logn_auto(), Etoile_bin::get_logn_comp(), Etoile::get_nbar(), Etoile::get_u_euler(), Cmp::integrale(), Etoile::is_relativistic(), mass_adm(), Etoile_bin::mass_b(), p_total_ener, star1, and star2.

double Binaire::virial (  )  const

Estimates the relative error on the virial theorem (for a relativistic one, it returns $|1 - M_{ Komar} / M_{ ADM}|$).

Definition at line 341 of file binaire_global.C.

References et, flat_scalar_prod(), Etoile::get_logn_auto(), Etoile_bin::get_logn_comp(), Etoile::get_nbar(), Etoile::get_press(), Etoile::get_u_euler(), Cmp::integrale(), Etoile::is_relativistic(), mass_adm(), mass_kom(), p_virial, star1, and star2.

double Binaire::virial_fus (  )  const

Estimates the relative error on the virial theorem calculated by J.L.Friedman, K.Uryu, and M.Shibata (PRD accepted, gr-qc/0108070).

The expression used in the LORENE is Eq.(5.7) in the paper by M.Shibata and K.Uryu (PRD64, 104017 (2001))

Definition at line 477 of file binaire_global.C.

References et, flat_scalar_prod(), Etoile::get_a_car(), Etoile_bin::get_akcar_auto(), Etoile_bin::get_akcar_comp(), Etoile_bin::get_d_beta_auto(), Etoile_bin::get_d_beta_comp(), Etoile_bin::get_d_logn_auto(), Etoile_bin::get_d_logn_comp(), Etoile::get_nnn(), Etoile::get_s_euler(), Cmp::integrale(), Etoile::is_relativistic(), mass_adm(), p_virial_fus, sqrt(), star1, star2, Cmp::std_base_scal(), and virial().

double Binaire::virial_gb (  )  const

Estimates the relative error on the virial theorem calculated by E.Gourgoulhon and S.Bonazzola (Class.

Quantum Grav. 11, 443 (1994): Eq.(29)) normalized by $2 G$.

Definition at line 404 of file binaire_global.C.

References et, flat_scalar_prod(), Etoile::get_a_car(), Etoile_bin::get_akcar_auto(), Etoile_bin::get_akcar_comp(), Etoile_bin::get_d_beta_auto(), Etoile_bin::get_d_beta_comp(), Etoile_bin::get_d_logn_auto(), Etoile_bin::get_d_logn_comp(), Etoile::get_s_euler(), Cmp::integrale(), Etoile::is_relativistic(), mass_adm(), p_virial_gb, sqrt(), star1, star2, Cmp::std_base_scal(), and virial().

void Binaire::write_global ( ostream &  ost  )  const

Friends And Related Function Documentation

ostream& operator<< ( ostream &  ,
const Binaire  
) [friend]

Save in a file.

Display


Member Data Documentation

Etoile_bin* Binaire::et[2] [private]

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 120 of file binaire.h.

double Binaire::omega [private]

Angular velocity with respect to an asymptotically inertial observer.

Definition at line 125 of file binaire.h.

Tbl* Binaire::p_angu_mom [mutable, private]

Total angular momentum of the system.

Definition at line 141 of file binaire.h.

double* Binaire::p_ham_constr [mutable, private]

Relative error on the Hamiltonian constraint.

Definition at line 156 of file binaire.h.

double* Binaire::p_mass_adm [mutable, private]

Total ADM mass of the system.

Definition at line 135 of file binaire.h.

double* Binaire::p_mass_kom [mutable, private]

Total Komar mass of the system.

Definition at line 138 of file binaire.h.

Tbl* Binaire::p_mom_constr [mutable, private]

Relative error on the momentum constraint.

Definition at line 159 of file binaire.h.

double* Binaire::p_total_ener [mutable, private]

Total energy of the system.

Definition at line 144 of file binaire.h.

double* Binaire::p_virial [mutable, private]

Virial theorem error.

Definition at line 147 of file binaire.h.

double* Binaire::p_virial_fus [mutable, private]

Virial theorem error by J.L.Friedman, K.Uryu, and M.Shibata.

Definition at line 153 of file binaire.h.

double* Binaire::p_virial_gb [mutable, private]

Virial theorem error by E.Gourgoulhon and S.Bonazzola.

Definition at line 150 of file binaire.h.

Cartesian triad of the absolute reference frame.

Definition at line 108 of file binaire.h.

First star of the system.

Definition at line 111 of file binaire.h.

Second star of the system.

Definition at line 114 of file binaire.h.

double Binaire::x_axe [private]

Absolute X coordinate of the rotation axis.

Definition at line 129 of file binaire.h.


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

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1