|
LORENE
|
Class for stars in binary system. More...
#include <star.h>
Public Member Functions | |
| Star_bin (Map &mp_i, int nzet_i, const Eos &eos_i, bool irrot, bool conf_flat) | |
| Standard constructor. | |
| Star_bin (const Star_bin &) | |
| Copy constructor. | |
| Star_bin (Map &mp_i, const Eos &eos_i, FILE *fich) | |
Constructor from a file (see sauve(FILE* )). | |
| virtual | ~Star_bin () |
| Destructor. | |
| void | operator= (const Star_bin &) |
Assignment to another Star_bin. | |
| Scalar & | set_pot_centri () |
| Read/write the centrifugal potential. | |
| Scalar & | set_logn_comp () |
| Read/write of the logarithm of the lapse generated principally by the companion. | |
| void | set_logn_auto (const Scalar &logn_auto_new) |
| Assignment of a new logn_auto. | |
| void | set_lnq_auto (const Scalar &lnq_auto_new) |
| Assignment of a new lnq_auto. | |
| Vector & | set_beta_auto () |
Read/write of ![]() | |
| Vector & | set_beta () |
Read/write of ![]() | |
| void | set_conf_flat (bool confflat) |
| Write if conformally flat. | |
| bool | is_irrotational () const |
Returns true for an irrotational motion, false for a corotating one. | |
| const Scalar & | get_psi0 () const |
| Returns the non-translational part of the velocity potential. | |
| const Vector & | get_d_psi () const |
| Returns the covariant derivative of the velocity potential (Spherical components with respect to the mapping of the star). | |
| const Vector & | get_wit_w () const |
| Returns the spatial projection of the fluid 3-velocity with respect to the co-orbiting observer. | |
| const Scalar & | get_loggam () const |
| Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer. | |
| const Vector & | get_bsn () const |
Returns the shift vector, divided by N, of the rotating coordinates, ![]() | |
| const Scalar & | get_pot_centri () const |
| Returns the centrifugal potential. | |
| const Scalar & | get_logn_auto () const |
| Returns the part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by the star. | |
| const Scalar & | get_logn_comp () const |
| Returns the part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by the companion star. | |
| const Vector & | get_beta_auto () const |
Returns the part of the shift vector ![]() | |
| const Vector & | get_beta_comp () const |
Returns the part of the shift vector ![]() | |
| const Scalar & | get_lnq_auto () const |
Returns the part of the vector field ![]() | |
| const Scalar & | get_lnq_comp () const |
Returns the part of the vector field ![]() | |
| const Vector & | get_dcov_logn () const |
Returns the covariant derivative of ![]() | |
| const Vector & | get_dcon_logn () const |
Returns the contravariant derivative of ![]() | |
| const Vector & | get_dcov_phi () const |
Returns the covariant derivative of ![]() | |
| const Vector & | get_dcon_phi () const |
Returns the contravariant derivative of ![]() | |
| const Scalar & | get_psi4 () const |
Return the conformal factor ![]() | |
| const Metric & | get_flat () const |
| Return the flat metric defined on the mapping (Spherical components with respect to the mapping of the star). | |
| const Metric & | get_gtilde () const |
Return the conformal 3-metric ![]() | |
| const Sym_tensor & | get_hij () const |
Return the total deviation of the inverse conformal metric ![]() | |
| const Sym_tensor & | get_hij_auto () const |
Return the deviation of the inverse conformal metric ![]() | |
| const Sym_tensor & | get_hij_comp () const |
Return the deviation of the inverse conformal metric ![]() | |
| const Sym_tensor & | get_tkij_auto () const |
Returns the part of the extrinsic curvature tensor ![]() beta_auto. | |
| const Sym_tensor & | get_tkij_comp () const |
Returns the part of the extrinsic curvature tensor ![]() beta_comp. | |
| const Scalar & | get_kcar_auto () const |
Returns the part of ![]() beta_auto. | |
| const Scalar & | get_kcar_comp () const |
Returns the part of ![]() beta_comp. | |
| const Scalar | get_decouple () const |
Returns the function used to construct beta_auto from beta. | |
| bool | is_conf_flat () const |
| virtual void | sauve (FILE *) const |
| Save in a file. | |
| virtual double | mass_b () const |
| Baryon mass. | |
| virtual double | mass_g () const |
| Gravitational mass. | |
| virtual double | xa_barycenter () const |
| Absolute coordinate X of the barycenter of the baryon density,. | |
| virtual void | hydro_euler () |
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame, as well as wit_w and loggam. | |
| void | update_metric (const Star_bin &comp, double omega) |
| Computes metric coefficients from known potentials, when the companion is another star. | |
| void | update_metric (const Star_bin &comp, const Star_bin &star_prev, double relax, double omega) |
Same as update_metric(const Star_bin& ) but with relaxation. | |
| void | update_metric_der_comp (const Star_bin &comp, double omega) |
| Computes the derivative of metric functions related to the companion star. | |
| void | kinematics (double omega, double x_axe) |
Computes the quantities bsn and pot_centri. | |
| void | fait_d_psi () |
Computes the gradient of the total velocity potential ![]() | |
| void | extrinsic_curvature (double omega) |
Computes tkij_auto and akcar_auto from beta_auto, nn and Q. | |
| void | equilibrium (double ent_c, int mermax, int mermax_potvit, int mermax_poisson, double relax_poisson, double relax_potvit, double thres_adapt, Tbl &diff, double om) |
| Computes an equilibrium configuration. | |
| double | velocity_potential (int mermax, double precis, double relax) |
Computes the non-translational part of the velocity scalar potential ![]() | |
| void | relaxation (const Star_bin &star_prev, double relax_ent, double relax_met, int mer, int fmer_met) |
Performs a relaxation on ent, logn_auto, lnq_auto, beta_auto and hij_auto. | |
| void | test_K_Hi () const |
| Test if the gauge conditions we impose are well satisfied. | |
| void | helical (double omega) const |
| Test of the helical symmetry. | |
| Map & | set_mp () |
| Read/write of the mapping. | |
| void | set_enthalpy (const Scalar &) |
| Assignment of the enthalpy field. | |
| void | equation_of_state () |
| Computes the proper baryon and energy density, as well as pressure from the enthalpy. | |
| virtual void | equilibrium_spher (double ent_c, double precis=1.e-14, const Tbl *pent_limit=0x0) |
| Computes a spherical static configuration. | |
| const Map & | get_mp () const |
| Returns the mapping. | |
| int | get_nzet () const |
| Returns the number of domains occupied by the star. | |
| const Eos & | get_eos () const |
| Returns the equation of state. | |
| const Scalar & | get_ent () const |
| Returns the enthalpy field. | |
| const Scalar & | get_nbar () const |
| Returns the proper baryon density. | |
| const Scalar & | get_ener () const |
| Returns the proper total energy density. | |
| const Scalar & | get_press () const |
| Returns the fluid pressure. | |
| const Scalar & | get_ener_euler () const |
| Returns the total energy density with respect to the Eulerian observer. | |
| const Scalar & | get_s_euler () const |
| Returns the trace of the stress tensor in the Eulerian frame. | |
| const Scalar & | get_gam_euler () const |
| Returns the Lorentz factor between the fluid and Eulerian observers. | |
| const Vector & | get_u_euler () const |
| Returns the fluid 3-velocity with respect to the Eulerian observer. | |
| const Tensor & | get_stress_euler () const |
| Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer. | |
| const Scalar & | get_logn () const |
| Returns the logarithm of the lapse N. | |
| const Scalar & | get_nn () const |
| Returns the lapse function N. | |
| const Vector & | get_beta () const |
Returns the shift vector ![]() | |
| const Scalar & | get_lnq () const |
| const Metric & | get_gamma () const |
Returns the 3-metric ![]() | |
| double | ray_eq () const |
Coordinate radius at ![]() ![]() | |
| double | ray_eq_pis2 () const |
Coordinate radius at ![]() ![]() | |
| double | ray_eq_pi () const |
Coordinate radius at ![]() ![]() | |
| double | ray_eq_3pis2 () const |
Coordinate radius at ![]() ![]() | |
| double | ray_pole () const |
Coordinate radius at ![]() | |
| virtual const Itbl & | l_surf () const |
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on the surface at the collocation points in ![]() | |
| const Tbl & | xi_surf () const |
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate ![]() ![]() | |
Protected Member Functions | |
| virtual void | del_deriv () const |
| Deletes all the derived quantities. | |
| virtual void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. | |
| virtual void | del_hydro_euler () |
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer. | |
| virtual ostream & | operator>> (ostream &) const |
| Operator >> (virtual function called by the operator <<). | |
Protected Attributes | |
| bool | irrotational |
true for an irrotational star, false for a corotating one | |
| Scalar | psi0 |
Scalar potential ![]() | |
| Vector | d_psi |
Gradient of ![]() | |
| Vector | wit_w |
| Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer. | |
| Scalar | loggam |
| Logarithm of the Lorentz factor between the fluid and the co-orbiting observer. | |
| Vector | bsn |
3-vector shift, divided by N, of the rotating coordinates, ![]() | |
| Scalar | pot_centri |
| Centrifugal potential. | |
| Scalar | logn_auto |
| Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by the star. | |
| Scalar | logn_comp |
| Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by the companion star. | |
| Vector | dcov_logn |
| Covariant derivative of the total logarithm of the lapse. | |
| Vector | dcon_logn |
| Contravariant derivative of the total logarithm of the lapse. | |
| Scalar | lnq_auto |
Scalar field ![]() | |
| Scalar | lnq_comp |
Scalar field ![]() | |
| Scalar | psi4 |
Conformal factor ![]() | |
| Vector | dcov_phi |
| Covariant derivative of the logarithm of the conformal factor. | |
| Vector | dcon_phi |
| Contravariant derivative of the logarithm of the conformal factor. | |
| Metric_flat | flat |
| Flat metric defined on the mapping (Spherical components with respect to the mapping of the star) . | |
| Metric | gtilde |
Conformal metric ![]() | |
| Vector | beta_auto |
| Part of the shift vector generated principally by the star (Spherical components with respect to the mapping of the star). | |
| Vector | beta_comp |
| Part of the shift vector generated principally by the star (Spherical components with respect to the mapping of the star). | |
| Sym_tensor | hij |
Total deviation of the inverse conformal metric ![]() | |
| Sym_tensor | hij_auto |
Deviation of the inverse conformal metric ![]() | |
| Sym_tensor | hij_comp |
Deviation of the inverse conformal metric ![]() | |
| Sym_tensor | tkij_auto |
Part of the extrinsic curvature tensor ![]() beta_auto. | |
| Sym_tensor | tkij_comp |
Part of the extrinsic curvature tensor ![]() beta_comp. | |
| Scalar | kcar_auto |
Part of the scalar ![]() beta_auto, i.e. | |
| Scalar | kcar_comp |
Part of the scalar ![]() beta_auto and beta_comp, i.e. | |
| Scalar | ssjm1_logn |
Effective source at the previous step for the resolution of the Poisson equation for logn_auto. | |
| Scalar | ssjm1_lnq |
Effective source at the previous step for the resolution of the Poisson equation for lnq_auto. | |
| Scalar | ssjm1_khi |
Effective source at the previous step for the resolution of the Poisson equation for khi. | |
| Vector | ssjm1_wbeta |
| Scalar | ssjm1_h11 |
Effective source at the previous step for the resolution of the Poisson equation for h00_auto. | |
| Scalar | ssjm1_h21 |
Effective source at the previous step for the resolution of the Poisson equation for h10_auto. | |
| Scalar | ssjm1_h31 |
Effective source at the previous step for the resolution of the Poisson equation for h20_auto. | |
| Scalar | ssjm1_h22 |
Effective source at the previous step for the resolution of the Poisson equation for h11_auto. | |
| Scalar | ssjm1_h32 |
Effective source at the previous step for the resolution of the Poisson equation for h21_auto. | |
| Scalar | ssjm1_h33 |
Effective source at the previous step for the resolution of the Poisson equation for h22_auto. | |
| Scalar | decouple |
Function used to construct the part ![]() ![]() | |
| bool | conf_flat |
true if the 3-metric is conformally flat, false for a more general metric. | |
| double * | p_xa_barycenter |
| Absolute coordinate X of the barycenter of the baryon density. | |
| Map & | mp |
| Mapping associated with the star. | |
| int | nzet |
Number of domains of *mp occupied by the star. | |
| const Eos & | eos |
| Equation of state of the stellar matter. | |
| Scalar | ent |
| Log-enthalpy. | |
| Scalar | nbar |
| Baryon density in the fluid frame. | |
| Scalar | ener |
| Total energy density in the fluid frame. | |
| Scalar | press |
| Fluid pressure. | |
| Scalar | ener_euler |
| Total energy density in the Eulerian frame. | |
| Scalar | s_euler |
| Trace of the stress scalar in the Eulerian frame. | |
| Scalar | gam_euler |
| Lorentz factor between the fluid and Eulerian observers. | |
| Vector | u_euler |
| Fluid 3-velocity with respect to the Eulerian observer. | |
| Sym_tensor | stress_euler |
| Spatial part of the stress-energy tensor with respect to the Eulerian observer. | |
| Scalar | logn |
| Logarithm of the lapse N . | |
| Scalar | nn |
| Lapse function N . | |
| Vector | beta |
| Shift vector. | |
| Scalar | lnq |
| Metric | gamma |
| 3-metric | |
| double * | p_ray_eq |
Coordinate radius at ![]() ![]() | |
| double * | p_ray_eq_pis2 |
Coordinate radius at ![]() ![]() | |
| double * | p_ray_eq_pi |
Coordinate radius at ![]() ![]() | |
| double * | p_ray_eq_3pis2 |
Coordinate radius at ![]() ![]() | |
| double * | p_ray_pole |
Coordinate radius at ![]() | |
| Itbl * | p_l_surf |
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in ![]() | |
| Tbl * | p_xi_surf |
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate ![]() ![]() | |
| double * | p_mass_b |
| Baryon mass. | |
| double * | p_mass_g |
| Gravitational mass. | |
Friends | |
| class | Binary |
Class for stars in binary system.
*** UNDER DEEVELOPMENT ***()
A Star_bin can be construted in two states, represented by the bool member irrotational: (i) irrotational (i.e. the fluid motion is irrotational) or (ii) rigidly corotating with respect to the orbital motion (synchronized binary).
| Lorene::Star_bin::Star_bin | ( | Map & | mp_i, |
| int | nzet_i, | ||
| const Eos & | eos_i, | ||
| bool | irrot, | ||
| bool | conf_flat ) |
Standard constructor.
| mp_i | Mapping on which the star will be defined |
| nzet_i | Number of domains occupied by the star |
| eos_i | Equation of state of the stellar matter |
| irrot | should be true for an irrotational star, false for a corotating one |
| conf_flat | should be true for a conformally flat metric false for a general one |
Definition at line 115 of file star_bin.C.
References Lorene::Star::beta, beta_auto, beta_comp, bsn, conf_flat, d_psi, dcon_logn, dcon_phi, dcov_logn, dcov_phi, decouple, flat, Lorene::Star::gamma, Lorene::get_bvect_cart(), gtilde, hij, hij_auto, hij_comp, irrotational, kcar_auto, kcar_comp, lnq_auto, lnq_comp, loggam, logn_auto, logn_comp, Lorene::Map(), Lorene::Star::mp, pot_centri, psi0, psi4, set_der_0x0(), ssjm1_h11, ssjm1_h21, ssjm1_h22, ssjm1_h31, ssjm1_h32, ssjm1_h33, ssjm1_khi, ssjm1_lnq, ssjm1_logn, Lorene::Star::Star(), Lorene::Star::stress_euler, tkij_auto, tkij_comp, Lorene::Star::u_euler, and wit_w.
| Lorene::Star_bin::Star_bin | ( | const Star_bin & | star | ) |
Copy constructor.
Definition at line 210 of file star_bin.C.
References beta_auto, beta_comp, bsn, conf_flat, d_psi, dcon_logn, dcon_phi, dcov_logn, dcov_phi, decouple, flat, gtilde, hij, hij_auto, hij_comp, irrotational, kcar_auto, kcar_comp, lnq_auto, lnq_comp, loggam, logn_auto, logn_comp, pot_centri, psi0, psi4, set_der_0x0(), ssjm1_h11, ssjm1_h21, ssjm1_h22, ssjm1_h31, ssjm1_h32, ssjm1_h33, ssjm1_khi, ssjm1_lnq, ssjm1_logn, Lorene::Star::Star(), Star_bin(), tkij_auto, tkij_comp, and wit_w.
Constructor from a file (see sauve(FILE* )).
| mp_i | Mapping on which the star will be defined |
| eos_i | Equation of state of the stellar matter |
| fich | input file (must have been created by the function sauve) |
Definition at line 258 of file star_bin.C.
References Lorene::Star::beta, beta_auto, beta_comp, bsn, conf_flat, d_psi, dcon_logn, dcon_phi, dcov_logn, dcov_phi, decouple, flat, Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::get_bvect_cart(), Lorene::get_mg(), gtilde, hij, hij_auto, hij_comp, irrotational, kcar_auto, kcar_comp, lnq_auto, lnq_comp, loggam, logn_auto, logn_comp, Lorene::Map(), Lorene::Star::mp, pot_centri, psi0, psi4, set_der_0x0(), ssjm1_h11, ssjm1_h21, ssjm1_h22, ssjm1_h31, ssjm1_h32, ssjm1_h33, ssjm1_khi, ssjm1_lnq, ssjm1_logn, Lorene::Star::Star(), Lorene::Star::stress_euler, tkij_auto, tkij_comp, Lorene::Star::u_euler, and wit_w.
|
virtual |
|
protectedvirtual |
Deletes all the derived quantities.
Reimplemented from Lorene::Star.
Definition at line 369 of file star_bin.C.
References Lorene::Star::del_deriv(), p_xa_barycenter, and set_der_0x0().
|
protectedvirtual |
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Reimplemented from Lorene::Star.
Definition at line 389 of file star_bin.C.
References del_deriv(), and Lorene::Star::del_hydro_euler().
|
inherited |
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition at line 462 of file star.C.
References Lorene::Param::add_int(), Lorene::Scalar::allocate_all(), del_deriv(), ener, ent, eos, Lorene::Mg3d::get_grille3d(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), mp, nbar, nzet, press, Lorene::Mtbl::set(), Lorene::Scalar::set_domain(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::std_spectral_base(), Lorene::Mtbl::t, and Lorene::Grille3d::x.
| void Lorene::Star_bin::equilibrium | ( | double | ent_c, |
| int | mermax, | ||
| int | mermax_potvit, | ||
| int | mermax_poisson, | ||
| double | relax_poisson, | ||
| double | relax_potvit, | ||
| double | thres_adapt, | ||
| Tbl & | diff, | ||
| double | om ) |
Computes an equilibrium configuration.
| ent_c | [input] Central enthalpy |
| mermax | [input] Maximum number of steps |
| mermax_poisson | [input] Maximum number of steps in poisson scalar |
| relax_poisson | [input] Relaxation factor in poisson scalar |
| mermax_potvit | [input] Maximum number of steps in Map_radial::poisson_compact |
| relax_potvit | [input] Relaxation factor in Map_radial::poisson_compact |
| thres_adapt | [input] Threshold on dH/dr for the adaptation of the mapping |
| diff | [output] 1-D Tbl for the storage of some error indicators |
Definition at line 143 of file star_bin_equilibrium.C.
References Lorene::abs(), Lorene::Param::add_cmp_mod(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Param::add_tbl(), Lorene::Param::add_tenseur_mod(), Lorene::Cmp::annule_hard(), Lorene::Scalar::annule_hard(), Lorene::Star::beta, beta_auto, conf_flat, Lorene::contract(), dcon_logn, dcon_phi, dcov_logn, dcov_phi, Lorene::Scalar::dec_dzpuis(), Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Tensor::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::diffrel(), Lorene::Vector::divergence(), Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::Star::equation_of_state(), Lorene::exp(), Lorene::Scalar::filtre(), flat, Lorene::Connection::get_delta(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), gtilde, hij, hij_auto, Lorene::Map_et::homothetie(), hydro_euler(), Lorene::Tensor::inc_dzpuis(), irrotational, kcar_auto, kcar_comp, lnq_auto, loggam, Lorene::Star::logn, logn_auto, logn_comp, Lorene::max(), Lorene::Star::mp, Lorene::Star::nn, Lorene::norme(), Lorene::Star::nzet, Lorene::phi, Lorene::Scalar::poisson(), Lorene::Tenseur::poisson_vect_oohara(), pot_centri, Lorene::pow(), Lorene::Star::press, psi4, Lorene::Star::ray_eq(), Lorene::Star::ray_eq_pi(), Lorene::Star::ray_pole(), Lorene::Star::s_euler, Lorene::Itbl::set(), Lorene::Tbl::set(), Lorene::Tenseur::set(), Lorene::Tensor::set(), Lorene::Vector::set(), Lorene::Cmp::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::set_outer_boundary(), Lorene::Tenseur::set_std_base(), Lorene::sqrt(), ssjm1_h11, ssjm1_h21, ssjm1_h22, ssjm1_h31, ssjm1_h32, ssjm1_h33, ssjm1_khi, ssjm1_lnq, ssjm1_logn, Lorene::Scalar::std_spectral_base(), Lorene::Tensor::std_spectral_base(), Lorene::Vector::std_spectral_base(), Lorene::Star::stress_euler, tkij_auto, tkij_comp, Lorene::Star::u_euler, Lorene::Scalar::val_grid_point(), and velocity_potential().
|
virtualinherited |
Computes a spherical static configuration.
| ent_c | [input] central value of the enthalpy |
| precis | [input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14) |
| ent_limit | [input] : array of enthalpy values to be set at the boundaries between the domains; if set to 0x0 (default), the initial values will be kept. |
Definition at line 95 of file star_equil_spher.C.
References Lorene::Map_et::adapt(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Param::add_tbl(), Lorene::diffrel(), Lorene::Map_af::dsdr(), Lorene::Scalar::dsdr(), ener, ener_euler, ent, equation_of_state(), Lorene::exp(), gam_euler, gamma, Lorene::Map_af::get_alpha(), Lorene::Map_et::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Map_et::get_beta(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map_af::homothetie(), Lorene::Scalar::integrale(), logn, mass_b(), mass_g(), mp, nn, Lorene::norme(), nzet, Lorene::Map_af::poisson(), press, s_euler, Lorene::Tensor::set(), Lorene::Map_af::set_alpha(), Lorene::Map_af::set_beta(), Lorene::Scalar::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Tensor::set_etat_zero(), Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), u_euler, and Lorene::Scalar::val_grid_point().
| void Lorene::Star_bin::extrinsic_curvature | ( | double | omega | ) |
Computes tkij_auto and akcar_auto from beta_auto, nn and Q.
| omega | angular velocity with respect to an asymptotically inertial observer |
Definition at line 75 of file star_bin_extr_curv.C.
References Lorene::Tensor::annule_domain(), Lorene::Scalar::annule_hard(), beta_auto, Lorene::contract(), del_deriv(), Lorene::Mg3d::get_nzone(), gtilde, hij_auto, kcar_auto, Lorene::Star::mp, Lorene::Star::nn, Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Scalar::set_spectral_va(), and tkij_auto.
| void Lorene::Star_bin::fait_d_psi | ( | ) |
Computes the gradient of the total velocity potential 
Definition at line 645 of file star_bin.C.
References bsn, Lorene::Vector::change_triad(), d_psi, Lorene::Star::ent, Lorene::exp(), flat, Lorene::Star::gam_euler, irrotational, Lorene::Star::mp, Lorene::Star::nzet, psi0, psi4, Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Vector::std_spectral_base(), and Lorene::Cmp::va.
|
inlineinherited |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inlineinherited |
|
inlineinherited |
Returns the total energy density with respect to the Eulerian observer.
Definition at line 376 of file star.h.
References ener_euler.
|
inlineinherited |
|
inlineinherited |
|
inline |
|
inlineinherited |
|
inlineinherited |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inlineinherited |
|
inline |
|
inline |
|
inline |
|
inlineinherited |
|
inline |
|
inline |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inline |
|
inlineinherited |
|
inline |
|
inline |
|
inlineinherited |
|
inlineinherited |
Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer.
Definition at line 390 of file star.h.
References stress_euler.
|
inline |
|
inline |
|
inlineinherited |
|
inline |
| void Lorene::Star_bin::helical | ( | double | omega | ) | const |
Test of the helical symmetry.
|
virtual |
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame, as well as wit_w and loggam.
The calculation is performed starting from the quantities ent, ener, press and bsn,
which are supposed to be up to date.
From these, the following fields are updated: gam_euler, u_euler, ener_euler, s_euler, stress_euler,
wit_w and loggam.
Reimplemented from Lorene::Star.
Definition at line 66 of file star_bin_hydro.C.
References bsn, Lorene::Tensor::change_triad(), Lorene::contract(), d_psi, del_deriv(), Lorene::Star::ener, Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::exp(), Lorene::Star::gam_euler, Lorene::Star::gamma, irrotational, Lorene::log(), loggam, Lorene::Star::mp, Lorene::norme(), Lorene::Star::press, Lorene::Star::s_euler, Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), Lorene::Star::stress_euler, Lorene::Star::u_euler, and wit_w.
|
inline |
Returns true for an irrotational motion, false for a corotating one.
Definition at line 773 of file star.h.
References irrotational.
| void Lorene::Star_bin::kinematics | ( | double | omega, |
| double | x_axe ) |
Computes the quantities bsn and pot_centri.
The calculation is performed starting from the quantities nn, beta, Q,
which are supposed to be up to date.
| omega | angular velocity with respect to an asymptotically inertial observer |
| x_axe | absolute X coordinate of the rotation axis |
Definition at line 75 of file star_bin_kinema.C.
References Lorene::Star::beta, bsn, Lorene::Tensor::change_triad(), Lorene::contract(), del_deriv(), Lorene::Star::gamma, Lorene::log(), Lorene::Star::mp, Lorene::Star::nn, pot_centri, Lorene::sqrt(), Lorene::xa, and Lorene::ya.
|
virtualinherited |
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on the surface at the collocation points in 
The stellar surface is defined as the location where the enthalpy (member ent) vanishes.
Reimplemented in Lorene::Star_rot.
Definition at line 63 of file star_global.C.
|
virtual |
Baryon mass.
Implements Lorene::Star.
Definition at line 67 of file star_bin_global.C.
References Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::Scalar::integrale(), Lorene::Star::nbar, Lorene::Star::p_mass_b, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
virtual |
Gravitational mass.
Implements Lorene::Star.
Definition at line 93 of file star_bin_global.C.
References Lorene::Star::ener_euler, Lorene::Star::gamma, Lorene::Scalar::integrale(), Lorene::Star::nn, Lorene::Star::p_mass_g, Lorene::Star::s_euler, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
| void Lorene::Star_bin::operator= | ( | const Star_bin & | star | ) |
Assignment to another Star_bin.
Definition at line 404 of file star_bin.C.
References beta_auto, beta_comp, bsn, conf_flat, d_psi, dcon_logn, dcon_phi, dcov_logn, dcov_phi, decouple, del_deriv(), flat, gtilde, hij, hij_auto, hij_comp, irrotational, kcar_auto, kcar_comp, lnq_auto, lnq_comp, loggam, logn_auto, logn_comp, Lorene::Star::operator=(), pot_centri, psi0, psi4, ssjm1_h11, ssjm1_h21, ssjm1_h22, ssjm1_h31, ssjm1_h32, ssjm1_h33, ssjm1_khi, ssjm1_lnq, ssjm1_logn, Star_bin(), tkij_auto, tkij_comp, and wit_w.
|
protectedvirtual |
Operator >> (virtual function called by the operator <<).
Reimplemented from Lorene::Star.
Definition at line 522 of file star_bin.C.
References Lorene::Star::beta, beta_auto, bsn, d_psi, Lorene::Star::gam_euler, irrotational, kcar_auto, kcar_comp, loggam, logn_auto, logn_comp, Lorene::max(), Lorene::min(), Lorene::Star::mp, Lorene::Star::operator>>(), Lorene::Star::ray_eq(), Lorene::Star::ray_eq_pi(), tkij_auto, tkij_comp, Lorene::Star::u_euler, wit_w, and xa_barycenter().
|
inherited |
Coordinate radius at 

Definition at line 108 of file star_global.C.
References Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), l_surf(), mp, p_ray_eq, Lorene::phi, and xi_surf().
|
inherited |
Coordinate radius at 

Definition at line 233 of file star_global.C.
References Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), l_surf(), mp, p_ray_eq_3pis2, Lorene::phi, ray_eq_pis2(), and xi_surf().
|
inherited |
Coordinate radius at 

Definition at line 186 of file star_global.C.
References Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), l_surf(), mp, p_ray_eq_pi, Lorene::phi, ray_eq(), and xi_surf().
|
inherited |
Coordinate radius at 

Definition at line 138 of file star_global.C.
References Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), l_surf(), mp, p_ray_eq_pis2, Lorene::phi, and xi_surf().
|
inherited |
Coordinate radius at 
Definition at line 278 of file star_global.C.
References Lorene::Mg3d::get_type_t(), l_surf(), mp, p_ray_pole, Lorene::phi, and xi_surf().
| void Lorene::Star_bin::relaxation | ( | const Star_bin & | star_prev, |
| double | relax_ent, | ||
| double | relax_met, | ||
| int | mer, | ||
| int | fmer_met ) |
Performs a relaxation on ent, logn_auto, lnq_auto, beta_auto and hij_auto.
| star_prev | [input] star at the previous step. |
| relax_ent | [input] Relaxation factor for ent |
| relax_met | [input] Relaxation factor for logn_auto, lnq_auto, beta_auto, only if (mer % fmer_met == 0). |
| mer | [input] Step number |
| fmer_met | [input] Step interval between metric updates |
Definition at line 701 of file star_bin.C.
References beta_auto, del_deriv(), Lorene::Star::ent, Lorene::Star::equation_of_state(), hij_auto, lnq_auto, logn_auto, and Star_bin().
|
virtual |
Save in a file.
Reimplemented from Lorene::Star.
Definition at line 489 of file star_bin.C.
References beta_auto, conf_flat, Lorene::Star::gam_euler, hij_auto, irrotational, lnq_auto, logn_auto, psi0, Lorene::Star::sauve(), ssjm1_h11, ssjm1_h21, ssjm1_h22, ssjm1_h31, ssjm1_h32, ssjm1_h33, ssjm1_khi, ssjm1_lnq, and ssjm1_logn.
| Vector & Lorene::Star_bin::set_beta | ( | ) |
Read/write of 
Definition at line 475 of file star_bin.C.
References Lorene::Star::beta, and del_deriv().
| Vector & Lorene::Star_bin::set_beta_auto | ( | ) |
|
inline |
|
protectedvirtual |
Sets to 0x0 all the pointers on derived quantities.
Reimplemented from Lorene::Star.
Definition at line 381 of file star_bin.C.
References p_xa_barycenter, and Lorene::Star::set_der_0x0().
|
inherited |
Assignment of the enthalpy field.
Definition at line 379 of file star.C.
References del_deriv(), ent, and equation_of_state().
|
inline |
|
inline |
| Scalar & Lorene::Star_bin::set_logn_comp | ( | ) |
Read/write of the logarithm of the lapse generated principally by the companion.
Definition at line 461 of file star_bin.C.
References del_deriv(), and logn_comp.
|
inlineinherited |
| Scalar & Lorene::Star_bin::set_pot_centri | ( | ) |
Read/write the centrifugal potential.
Definition at line 454 of file star_bin.C.
References del_deriv(), and pot_centri.
| void Lorene::Star_bin::test_K_Hi | ( | ) | const |
Test if the gauge conditions we impose are well satisfied.
Definition at line 727 of file star_bin.C.
References flat, gtilde, Lorene::Star::mp, and Lorene::norme().
| void Lorene::Star_bin::update_metric | ( | const Star_bin & | comp, |
| const Star_bin & | star_prev, | ||
| double | relax, | ||
| double | omega ) |
Same as update_metric(const Star_bin& ) but with relaxation.
| comp | companion star. |
| star_prev | previous value of the star. |
| relax | relaxation parameter. |
| omega | angular velocity with respect to an asymptotically inertial observer |
Definition at line 276 of file star_bin_upmetr.C.
References Lorene::Star::beta, beta_auto, beta_comp, Lorene::Tensor::change_triad(), Lorene::Vector::change_triad(), conf_flat, del_deriv(), Lorene::exp(), extrinsic_curvature(), flat, Lorene::Star::gamma, Lorene::Star::get_mp(), Lorene::Tensor::get_triad(), gtilde, hij, hij_auto, hij_comp, lnq_auto, lnq_comp, Lorene::Star::logn, logn_auto, logn_comp, Lorene::Map(), Lorene::Star::mp, Lorene::Star::nn, psi4, Lorene::Tensor::set(), Star_bin(), and Lorene::Scalar::val_grid_point().
| void Lorene::Star_bin::update_metric | ( | const Star_bin & | comp, |
| double | omega ) |
Computes metric coefficients from known potentials, when the companion is another star.
The calculation is performed starting from the quantities logn_auto, lnq_auto, beta_auto, hij_auto, comp.logn_auto, comp.lnq_auto, comp.beta_auto, comp.hij_auto which are supposed to be up to date. From these, the following fields are updated: logn_comp, lnq_comp, beta_comp, hij_comp, nn, psi4, beta,
| comp | companion star. |
| omega | angular velocity with respect to an asymptotically inertial observer |
Definition at line 93 of file star_bin_upmetr.C.
References Lorene::Star::beta, beta_auto, beta_comp, Lorene::Tensor::change_triad(), Lorene::Vector::change_triad(), conf_flat, del_deriv(), Lorene::exp(), extrinsic_curvature(), flat, Lorene::Star::gamma, Lorene::Star::get_mp(), Lorene::Tensor::get_triad(), gtilde, hij, hij_auto, hij_comp, lnq_auto, lnq_comp, Lorene::Star::logn, logn_auto, logn_comp, Lorene::Map(), Lorene::Star::mp, Lorene::Star::nn, psi4, Lorene::Tensor::set(), Star_bin(), and Lorene::Scalar::val_grid_point().
| void Lorene::Star_bin::update_metric_der_comp | ( | const Star_bin & | comp, |
| double | omega ) |
Computes the derivative of metric functions related to the companion star.
| omega | angular velocity with respect to an asymptotically inertial observer |
Definition at line 101 of file star_bin_upmetr_der.C.
References Lorene::Tensor::annule_domain(), Lorene::Scalar::annule_hard(), beta_comp, Lorene::Vector::change_triad(), Lorene::contract(), dcon_logn, dcon_phi, dcov_logn, dcov_phi, Lorene::Tensor::dec_dzpuis(), del_deriv(), flat, get_flat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_triad(), gtilde, hij_comp, kcar_comp, lnq_auto, logn_auto, Lorene::Star::mp, Lorene::Star::nn, Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Scalar::set_spectral_va(), Star_bin(), tkij_auto, and tkij_comp.
| double Lorene::Star_bin::velocity_potential | ( | int | mermax, |
| double | precis, | ||
| double | relax ) |
Computes the non-translational part of the velocity scalar potential 
| mermax | [input] Maximum number of steps in the iteration |
| precis | [input] Required precision: the iteration will be stopped when the relative difference on ![]() precis. |
| relax | [input] Relaxation factor. |
Definition at line 76 of file star_bin_vel_pot.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Scalar::annule(), bsn, Lorene::Vector::change_triad(), Lorene::contract(), d_psi, Lorene::Scalar::derive_con(), Lorene::Scalar::derive_cov(), Lorene::diffrel(), Lorene::Tensor::down(), Lorene::Star::ent, Lorene::Star::eos, Lorene::exp(), flat, Lorene::Star::gam_euler, Lorene::Tensor::get_triad(), hij, Lorene::Tensor::inc_dzpuis(), Lorene::Star::mp, Lorene::norme(), Lorene::Star::nzet, psi0, psi4, Lorene::Tenseur::set(), Lorene::Vector::set(), Lorene::Scalar::set_domain(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tensor::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Tensor::set_triad(), Lorene::Scalar::std_spectral_base(), Lorene::Vector::std_spectral_base(), Lorene::Cmp::va, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
virtual |
Absolute coordinate X of the barycenter of the baryon density,.
Definition at line 118 of file star_bin_global.C.
References Lorene::Tensor::annule_domain(), Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::Scalar::integrale(), mass_b(), Lorene::Star::mp, Lorene::Star::nbar, p_xa_barycenter, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
inherited |
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate 

The stellar surface is defined as the location where the enthalpy (member ent) vanishes.
Definition at line 89 of file star_global.C.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protected |
|
protectedinherited |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protectedinherited |
|
protected |
|
protected |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotected |
|
mutableprotectedinherited |
|
protected |
|
protectedinherited |
|
protected |
|
protectedinherited |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protectedinherited |
|
protected |
|
protected |
|
protectedinherited |
|
protected |