LORENE
Lorene::Binary Class Reference

Binary systems. More...

#include <binary.h>

Public Member Functions

 Binary (Map &mp1, int nzet1, const Eos &eos1, int irrot1, Map &mp2, int nzet2, const Eos &eos2, int irrot2, int conf_flat)
 Standard constructor. More...
 
 Binary (const Binary &)
 Copy constructor. More...
 
 Binary (Map &mp1, const Eos &eos1, Map &mp2, const Eos &eos2, FILE *fich)
 Constructor from a file (see sauve(FILE* )). More...
 
 ~Binary ()
 Destructor. More...
 
void operator= (const Binary &)
 Assignment to another Binary. More...
 
Star_binset (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_binoperator() (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 Tblangu_mom () const
 Total angular 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 ham_constr () const
 Estimates the relative error on the Hamiltonian constraint. More...
 
const Tblmom_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...
 
void fait_decouple ()
 Calculates decouple which is used to obtain qq_auto by the formula : qq_auto = decouple * qq. More...
 
void dirac_gauge ()
 Function used to impose Dirac gauge during an iteration. More...
 
void helical ()
 Function testing the helical symmetry. 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 star1
 First star of the system. More...
 
Star_bin star2
 Second star of the system. More...
 
Star_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. 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...
 
Tblp_angu_mom
 Total angular momentum of the system. More...
 
double * p_total_ener
 Total energy of the system. More...
 
double * p_virial
 Virial theorem error. More...
 
double * p_ham_constr
 Relative error on the Hamiltonian constraint. More...
 
Tblp_mom_constr
 Relative error on the momentum constraint. More...
 

Private Member Functions

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

Friends

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

Detailed Description

Binary systems.

*** UNDER DEVELOPMENT ***()

Version
#$Id: binary.h,v 1.10 2014/10/13 08:52:32 j_novak Exp $#

Definition at line 73 of file binary.h.

Constructor & Destructor Documentation

◆ Binary() [1/3]

Lorene::Binary::Binary ( Map mp1,
int  nzet1,
const Eos eos1,
int  irrot1,
Map mp2,
int  nzet2,
const Eos eos2,
int  irrot2,
int  conf_flat 
)

Standard constructor.

Parameters
mp1Mapping on which star1 will be defined
nzet1Number of domains occupied by star1
eos1Equation of state of star1
irrot1should be true if star1 is irrotational, false if star1 is corotating
mp2Mapping on which star2 will be defined
nzet2Number of domains occupied by star2
eos2Equation of state of star2
irrot2should be true if star2 is irrotational, false if star2 is corotating
conf_flatshould be true for a 3-metric conformally flat and false for a more general one.

Definition at line 102 of file binary.C.

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

◆ Binary() [2/3]

Lorene::Binary::Binary ( const Binary bibi)

Copy constructor.

Definition at line 121 of file binary.C.

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

◆ Binary() [3/3]

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

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

Parameters
mp1Mapping on which star1 will be defined
eos1Equation of state of star1
mp2Mapping on which star2 will be defined
eos2Equation of state of star2
fichinput file (must have been created by the function sauve)

Definition at line 136 of file binary.C.

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

◆ ~Binary()

Lorene::Binary::~Binary ( )

Destructor.

Definition at line 157 of file binary.C.

References del_deriv().

Member Function Documentation

◆ analytical_omega()

void Lorene::Binary::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 binary_omegaana.C.

References del_deriv(), Lorene::Star_bin::is_irrotational(), Lorene::Star_bin::mass_g(), omega, Lorene::pow(), Lorene::Star::ray_eq(), separation(), Lorene::sqrt(), star1, and star2.

◆ analytical_shift()

◆ angu_mom()

◆ del_deriv()

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

Deletes all the derived quantities.

Definition at line 167 of file binary.C.

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

◆ dirac_gauge()

◆ display_poly()

◆ fait_decouple()

void Lorene::Binary::fait_decouple ( )

◆ get_omega()

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

Returns the orbital angular velocity [f_unit].

Definition at line 207 of file binary.h.

References omega.

◆ get_x_axe()

double Lorene::Binary::get_x_axe ( ) const
inline

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

Definition at line 210 of file binary.h.

References x_axe.

◆ ham_constr()

double Lorene::Binary::ham_constr ( ) const

Estimates the relative error on the Hamiltonian constraint.

◆ helical()

void Lorene::Binary::helical ( )

Function testing the helical symmetry.

Definition at line 79 of file binary_helical.C.

References Lorene::Tbl::annule_hard(), Lorene::Star::beta, Lorene::Star_bin::beta_auto, Lorene::Star_bin::beta_comp, Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Star_bin::dcon_logn, Lorene::Star_bin::dcon_phi, Lorene::Star_bin::dcov_logn, Lorene::Star_bin::dcov_phi, Lorene::Scalar::dec_dzpuis(), Lorene::Tensor::derive_con(), Lorene::Scalar::derive_con(), Lorene::Tensor_sym::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Sym_tensor::derive_lie(), Lorene::Metric::determinant(), Lorene::Sym_tensor::divergence(), Lorene::Vector::divergence(), Lorene::Star::ener_euler, et, Lorene::exp(), Lorene::Star_bin::flat, Lorene::Star::gamma, Lorene::Map::get_bvect_cart(), Lorene::Star::get_gamma(), Lorene::Map::get_mg(), Lorene::Star::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Map::get_rot_phi(), Lorene::Star_bin::gtilde, Lorene::Star_bin::hij, Lorene::Star_bin::hij_auto, Lorene::Star_bin::hij_comp, Lorene::Tenseur::inc_dzpuis(), Lorene::Tensor::inc_dzpuis(), Lorene::Scalar::inc_dzpuis(), Lorene::Scalar::integrale_domains(), Lorene::Scalar::laplacian(), Lorene::Star_bin::lnq_auto, Lorene::Star::logn, Lorene::Star_bin::logn_auto, mass_adm(), Lorene::Star::mp, Lorene::Star::nn, omega, Lorene::pow(), Lorene::Star_bin::psi4, Lorene::Map::r, Lorene::Star::s_euler, Lorene::Coord::set(), Lorene::Itbl::set(), Lorene::Tbl::set(), Lorene::Tensor::set(), Lorene::sqrt(), star1, star2, Lorene::Scalar::std_spectral_base(), Lorene::Star::stress_euler, Lorene::Star_bin::tkij_auto, Lorene::Star_bin::tkij_comp, Lorene::Tensor::up_down(), Lorene::Map::val_r(), Lorene::Map::xa, Lorene::Map::ya, and Lorene::Map::za.

◆ mass_adm()

◆ mass_adm_vol()

◆ mass_kom()

◆ mass_kom_vol()

◆ mom_constr()

const Tbl& Lorene::Binary::mom_constr ( ) const

Estimates the relative error on the momentum constraint.

◆ operator()()

const Star_bin& Lorene::Binary::operator() ( int  i) const
inline

Returns a reference to the star no. i.

Definition at line 202 of file binary.h.

References et.

◆ operator=()

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

Assignment to another Binary.

Definition at line 203 of file binary.C.

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

◆ operator>>()

ostream & Lorene::Binary::operator>> ( ostream &  ost) const
private

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

Definition at line 239 of file binary.C.

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

◆ orbit()

void Lorene::Binary::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.

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)
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 90 of file binary_orbite.C.

References Lorene::Param::add_double(), Lorene::Vector::change_triad(), Lorene::Tensor::change_triad(), Lorene::Metric::cov(), Lorene::Scalar::dsdx(), et, Lorene::Star::get_beta(), Lorene::Star::get_ent(), Lorene::Star::get_eos(), Lorene::Star::get_gamma(), Lorene::Star_bin::get_loggam(), Lorene::Star_bin::get_logn_auto(), Lorene::Star_bin::get_logn_comp(), Lorene::Star::get_mp(), Lorene::Star::get_nn(), Lorene::Tbl::get_taille(), Lorene::norme(), omega, Lorene::Scalar::val_grid_point(), x_axe, Lorene::Star_bin::xa_barycenter(), Lorene::zero_list(), Lorene::zerosec(), and Lorene::zerosec_b().

◆ sauve()

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

Save in a file.

Definition at line 221 of file binary.C.

References Lorene::fwrite_be(), omega, Lorene::Star_bin::sauve(), star1, star2, and x_axe.

◆ separation()

double Lorene::Binary::separation ( ) const

Returns the coordinate separation of the two stellar centers [r_unit].

Definition at line 610 of file binary.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.

◆ set()

Star_bin& Lorene::Binary::set ( int  i)
inline

Read/write of the star no. i.

Definition at line 186 of file binary.h.

References del_deriv(), and et.

◆ set_der_0x0()

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

Sets to 0x0 all the pointers on derived quantities.

Definition at line 183 of file binary.C.

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

◆ set_omega()

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

Sets the orbital angular velocity [f_unit].

Definition at line 192 of file binary.h.

References omega.

◆ set_x_axe()

double& Lorene::Binary::set_x_axe ( )
inline

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

Definition at line 195 of file binary.h.

References x_axe.

◆ total_ener()

double Lorene::Binary::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_{\rm ADM} - M_{\rm bar,1} - M_{\rm bar,2}$.

Definition at line 610 of file binary_global.C.

References p_total_ener.

◆ virial()

double Lorene::Binary::virial ( ) const

Estimates the relative error on the virial theorem.

Definition at line 630 of file binary_global.C.

References mass_adm(), mass_kom(), and p_virial.

◆ write_global()

Friends And Related Function Documentation

◆ operator<<

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

Display.

Definition at line 233 of file binary.C.

Member Data Documentation

◆ et

Star_bin* Lorene::Binary::et[2]
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 89 of file binary.h.

◆ omega

double Lorene::Binary::omega
protected

Angular velocity with respect to an asymptotically inertial observer.

Definition at line 94 of file binary.h.

◆ p_angu_mom

Tbl* Lorene::Binary::p_angu_mom
mutableprotected

Total angular momentum of the system.

Definition at line 111 of file binary.h.

◆ p_ham_constr

double* Lorene::Binary::p_ham_constr
mutableprotected

Relative error on the Hamiltonian constraint.

Definition at line 120 of file binary.h.

◆ p_mass_adm

double* Lorene::Binary::p_mass_adm
mutableprotected

Total ADM mass of the system.

Definition at line 105 of file binary.h.

◆ p_mass_kom

double* Lorene::Binary::p_mass_kom
mutableprotected

Total Komar mass of the system.

Definition at line 108 of file binary.h.

◆ p_mom_constr

Tbl* Lorene::Binary::p_mom_constr
mutableprotected

Relative error on the momentum constraint.

Definition at line 123 of file binary.h.

◆ p_total_ener

double* Lorene::Binary::p_total_ener
mutableprotected

Total energy of the system.

Definition at line 114 of file binary.h.

◆ p_virial

double* Lorene::Binary::p_virial
mutableprotected

Virial theorem error.

Definition at line 117 of file binary.h.

◆ star1

Star_bin Lorene::Binary::star1
protected

First star of the system.

Definition at line 80 of file binary.h.

◆ star2

Star_bin Lorene::Binary::star2
protected

Second star of the system.

Definition at line 83 of file binary.h.

◆ x_axe

double Lorene::Binary::x_axe
protected

Absolute X coordinate of the rotation axis.

Definition at line 98 of file binary.h.


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