|
LORENE
|
Class for relativistic differentially rotating stars in Dirac gauge and maximal slicing. More...
#include <star_rot_dirac_diff.h>
Public Member Functions | |
| Star_rot_Dirac_diff (Map &mp_i, int nzet_i, const Eos &eos_i, double(*frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i) | |
| Standard constructor. More... | |
| Star_rot_Dirac_diff (const Star_rot_Dirac_diff &) | |
| Copy constructor. More... | |
| Star_rot_Dirac_diff (Map &mp_i, const Eos &eos_i, FILE *fich, double(*frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &)) | |
Constructor from a file (see sauve(FILE*) ). More... | |
| virtual | ~Star_rot_Dirac_diff () |
| Destructor. More... | |
| void | operator= (const Star_rot_Dirac_diff &) |
Assignment to another Star_rot_Dirac_diff. More... | |
| const Scalar & | get_omega_field () const |
Returns the angular velocity field . More... | |
| virtual double | get_omega_c () const |
Returns the central value of the rotation angular velocity ([f_unit] ) More... | |
| virtual void | sauve (FILE *) const |
| Save in a file. More... | |
| virtual double | tsw () const |
| Ratio T/W. More... | |
| virtual void | hydro_euler () |
| Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame. More... | |
| void | fait_omega_field (double omeg_min, double omeg_max, double precis, int nitermax) |
Computes (member omega_field ). More... | |
| void | fait_prim_field () |
Computes the member prim_field from omega_field . More... | |
| double | funct_omega (double omeg) const |
Evaluates , where F is the function defining the rotation profile. More... | |
| double | prim_funct_omega (double omeg) const |
Evaluates the primitive of , where F is the function defining the rotation profile. More... | |
| virtual void | equilibrium (double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff) |
| Computes an equilibrium configuration. More... | |
| int | get_relat () const |
| Returns the relativity parameter. More... | |
| bool | is_relativistic () const |
| Checks whether the star is computed using a relativistic theory. More... | |
| int | spectral_filter_order () const |
| Returns the filtering order. More... | |
| double | get_omega () const |
Returns the rotation angular velocity . More... | |
| const Scalar & | get_psi4 () const |
Returns the conformal factor . More... | |
| const Scalar & | get_psi2 () const |
Returns . More... | |
| const Scalar & | get_qqq () const |
Returns . More... | |
| const Scalar & | get_ln_psi () const |
Returns . More... | |
| const Vector & | get_j_euler () const |
| Returns the momentum density 3-vector with respect to the Eulerian observer. More... | |
| const Scalar & | get_v2 () const |
Returns . More... | |
| const Metric | get_tgamma () const |
Returns the conformal metric . More... | |
| const Sym_tensor | get_aa () const |
Returns . More... | |
| const Sym_tensor | get_taa () const |
Returns . More... | |
| const Scalar | get_aa_quad () const |
Returns . More... | |
| const Sym_tensor_trans | get_hh () const |
Returns . More... | |
| virtual double | mass_b () const |
| Baryonic mass. More... | |
| virtual double | mass_g () const |
| Gravitational mass. More... | |
| virtual double | angu_mom () const |
| Angular momentum. More... | |
| virtual double | grv2 () const |
| Error on the virial identity GRV2. More... | |
| virtual double | grv3 () const |
| Error on the virial identity GRV3. More... | |
| virtual double | aplat () const |
| Flattening r_pole/r_eq. More... | |
| virtual double | r_circ () const |
| Circumferential equatorial radius. More... | |
| virtual double | rp_circ () const |
| Circumferential polar radius. More... | |
| virtual double | ellipt () const |
| Ellipticity e. More... | |
| void | update_metric () |
| Computes metric quantities from known potentials. More... | |
| void | solve_logn_f (Scalar &ln_f_new) const |
| Solution of the two scalar Poisson equations for rotating stars in Dirac gauge. More... | |
| void | solve_logn_q (Scalar &ln_q_new) const |
| Solution of the two scalar Poisson equations for rotating stars in Dirac gauge. More... | |
| void | solve_qqq (Scalar &q_new) const |
| Solution of the two scalar Poisson equations for rotating stars in Dirac gauge. More... | |
| void | solve_shift (Vector &shift_new) const |
| Solution of the shift equation for rotating stars in Dirac gauge. More... | |
| void | solve_hij (Sym_tensor_trans &hij_new) const |
| Solution of the tensor Poisson equation for rotating stars in Dirac gauge. More... | |
| Map & | set_mp () |
| Read/write of the mapping. More... | |
| void | set_enthalpy (const Scalar &) |
| Assignment of the enthalpy field. More... | |
| void | equation_of_state () |
| Computes the proper baryon and energy density, as well as pressure from the enthalpy. More... | |
| virtual void | equilibrium_spher (double ent_c, double precis=1.e-14, const Tbl *pent_limit=0x0) |
| Computes a spherical static configuration. More... | |
| const Map & | get_mp () const |
| Returns the mapping. More... | |
| int | get_nzet () const |
| Returns the number of domains occupied by the star. More... | |
| const Eos & | get_eos () const |
| Returns the equation of state. More... | |
| const Scalar & | get_ent () const |
| Returns the enthalpy field. More... | |
| const Scalar & | get_nbar () const |
| Returns the proper baryon density. More... | |
| const Scalar & | get_ener () const |
| Returns the proper total energy density. More... | |
| const Scalar & | get_press () const |
| Returns the fluid pressure. More... | |
| const Scalar & | get_ener_euler () const |
| Returns the total energy density with respect to the Eulerian observer. More... | |
| const Scalar & | get_s_euler () const |
| Returns the trace of the stress tensor in the Eulerian frame. More... | |
| const Scalar & | get_gam_euler () const |
| Returns the Lorentz factor between the fluid and Eulerian observers. More... | |
| const Vector & | get_u_euler () const |
| Returns the fluid 3-velocity with respect to the Eulerian observer. More... | |
| const Tensor & | get_stress_euler () const |
| Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer. More... | |
| const Scalar & | get_logn () const |
| Returns the logarithm of the lapse N. More... | |
| const Scalar & | get_nn () const |
| Returns the lapse function N. More... | |
| const Vector & | get_beta () const |
Returns the shift vector . More... | |
| const Scalar & | get_lnq () const |
| const Metric & | get_gamma () const |
Returns the 3-metric . More... | |
| double | ray_eq () const |
Coordinate radius at , [r_unit]. More... | |
| double | ray_eq_pis2 () const |
Coordinate radius at , [r_unit]. More... | |
| double | ray_eq_pi () const |
Coordinate radius at , [r_unit]. More... | |
| double | ray_eq_3pis2 () const |
Coordinate radius at , [r_unit]. More... | |
| double | ray_pole () const |
Coordinate radius at [r_unit]. More... | |
| 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 . More... | |
| const Tbl & | xi_surf () const |
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate on the surface at the collocation points in . More... | |
Protected Member Functions | |
| virtual ostream & | operator>> (ostream &) const |
| Operator >> (virtual function called by the operator <<). More... | |
| virtual 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... | |
| virtual void | del_hydro_euler () |
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer. More... | |
Protected Attributes | |
| double(* | frot )(double, const Tbl &) |
Function defining the rotation profile. More... | |
| double(* | primfrot )(double, const Tbl &) |
Primitive of the function , which vanishes at the stellar center. More... | |
| Tbl | par_frot |
Parameters of the function . More... | |
| Scalar | omega_field |
Field . More... | |
| double | omega_min |
Minimum value of . More... | |
| double | omega_max |
Maximum value of . More... | |
| Scalar | prim_field |
Field . More... | |
| int | relat_type |
| Relativistic flag. More... | |
| int | spectral_filter |
| Spectral exponential filtering order. More... | |
| double | omega |
Rotation angular velocity ([f_unit] ) More... | |
| Scalar | psi4 |
Conformal factor . More... | |
| Scalar | psi2 |
More... | |
| Scalar | qqq |
More... | |
| Scalar | ln_psi |
More... | |
| Vector | j_euler |
| Momentum density 3-vector with respect to the Eulerian observer. More... | |
| Scalar | v2 |
More... | |
| const Metric_flat & | flat |
flat metric (spherical components) More... | |
| Metric | tgamma |
More... | |
| Sym_tensor | aa |
More... | |
| Sym_tensor | taa |
More... | |
| Scalar | aa_quad |
More... | |
| Sym_tensor_trans | hh |
is defined by . More... | |
| double * | p_angu_mom |
| Angular momentum. More... | |
| double * | p_grv2 |
| Error on the virial identity GRV2. More... | |
| double * | p_grv3 |
| Error on the virial identity GRV3. More... | |
| double * | p_tsw |
| Ratio T/W. More... | |
| double * | p_r_circ |
| Circumferential equatorial radius. More... | |
| double * | p_rp_circ |
| Circumferential polar radius. More... | |
| Map & | mp |
| Mapping associated with the star. More... | |
| int | nzet |
Number of domains of *mp occupied by the star. More... | |
| const Eos & | eos |
| Equation of state of the stellar matter. More... | |
| Scalar | ent |
| Log-enthalpy. More... | |
| Scalar | nbar |
| Baryon density in the fluid frame. More... | |
| Scalar | ener |
| Total energy density in the fluid frame. More... | |
| Scalar | press |
| Fluid pressure. More... | |
| Scalar | ener_euler |
| Total energy density in the Eulerian frame. More... | |
| Scalar | s_euler |
| Trace of the stress scalar in the Eulerian frame. More... | |
| Scalar | gam_euler |
| Lorentz factor between the fluid and Eulerian observers. More... | |
| Vector | u_euler |
| Fluid 3-velocity with respect to the Eulerian observer. More... | |
| Sym_tensor | stress_euler |
| Spatial part of the stress-energy tensor with respect to the Eulerian observer. More... | |
| Scalar | logn |
| Logarithm of the lapse N . More... | |
| Scalar | nn |
| Lapse function N . More... | |
| Vector | beta |
| Shift vector. More... | |
| Scalar | lnq |
| Metric | gamma |
| 3-metric More... | |
| double * | p_ray_eq |
Coordinate radius at , . More... | |
| double * | p_ray_eq_pis2 |
Coordinate radius at , . More... | |
| double * | p_ray_eq_pi |
Coordinate radius at , . More... | |
| double * | p_ray_eq_3pis2 |
Coordinate radius at , . More... | |
| double * | p_ray_pole |
Coordinate radius at . More... | |
| 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 . More... | |
| Tbl * | p_xi_surf |
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the surface at the collocation points in . More... | |
| double * | p_mass_b |
| Baryon mass. More... | |
| double * | p_mass_g |
| Gravitational mass. More... | |
Friends | |
| ostream & | operator<< (ostream &, const Star &) |
| Display. More... | |
Class for relativistic differentially rotating stars in Dirac gauge and maximal slicing.
(*** Under development ***)()
Definition at line 49 of file star_rot_dirac_diff.h.
| Lorene::Star_rot_Dirac_diff::Star_rot_Dirac_diff | ( | Map & | mp_i, |
| int | nzet_i, | ||
| const Eos & | eos_i, | ||
| double(*)(double, const Tbl &) | frot_i, | ||
| double(*)(double, const Tbl &) | primfrot_i, | ||
| const Tbl & | par_frot_i | ||
| ) |
Standard constructor.
| mp_i | Mapping on which the star will be defined |
| nzet_i | Number of domains occupied by the star |
| relat | should be true for a relativistic star, false for a Newtonian one |
| eos_i | Equation of state of the stellar matter |
| frot_i | Function defining the rotation profile. |
| primfrot_i | Primitive of which vanishes at the stellar center |
| par_frot_i | Parameters of functions frot_i and primfrot_i , par_frot_i(0) being the central value of . |
Definition at line 53 of file star_rot_dirac_diff.C.
References omega_field, omega_max, omega_min, and prim_field.
| Lorene::Star_rot_Dirac_diff::Star_rot_Dirac_diff | ( | const Star_rot_Dirac_diff & | star | ) |
Copy constructor.
Definition at line 76 of file star_rot_dirac_diff.C.
| Lorene::Star_rot_Dirac_diff::Star_rot_Dirac_diff | ( | Map & | mp_i, |
| const Eos & | eos_i, | ||
| FILE * | fich, | ||
| double(*)(double, const Tbl &) | frot_i, | ||
| double(*)(double, const Tbl &) | primfrot_i | ||
| ) |
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 ) |
| frot_i | Function defining the rotation profile. |
| primfrot_i | Primitive of which vanishes at the stellar center |
Definition at line 90 of file star_rot_dirac_diff.C.
References fait_prim_field(), Lorene::fread_be(), Lorene::Map::get_mg(), Lorene::Star::mp, omega_field, omega_max, and omega_min.
|
virtual |
Destructor.
Definition at line 116 of file star_rot_dirac_diff.C.
|
virtualinherited |
Angular momentum.
Definition at line 146 of file strot_dirac_global.C.
References Lorene::contract(), Lorene::Metric::cov(), Lorene::Metric::determinant(), Lorene::Star::gamma, Lorene::Map::get_bvect_spher(), Lorene::Scalar::integrale(), Lorene::Star_rot_Dirac::j_euler, Lorene::Star::mp, Lorene::Star_rot_Dirac::p_angu_mom, Lorene::Vector::set(), Lorene::Scalar::set_etat_zero(), Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
virtualinherited |
Flattening r_pole/r_eq.
Definition at line 421 of file strot_dirac_global.C.
References Lorene::Star::ray_eq(), and Lorene::Star::ray_pole().
|
protectedvirtualinherited |
Deletes all the derived quantities.
Reimplemented from Lorene::Star.
Definition at line 235 of file star_rot_dirac.C.
References Lorene::Star::del_deriv(), Lorene::Star_rot_Dirac::p_angu_mom, Lorene::Star_rot_Dirac::p_grv2, Lorene::Star_rot_Dirac::p_grv3, Lorene::Star_rot_Dirac::p_r_circ, Lorene::Star_rot_Dirac::p_rp_circ, Lorene::Star_rot_Dirac::p_tsw, and Lorene::Star_rot_Dirac::set_der_0x0().
|
protectedvirtualinherited |
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Reimplemented from Lorene::Star.
Definition at line 263 of file star_rot_dirac.C.
References Lorene::Star_rot_Dirac::del_deriv(), Lorene::Star::del_hydro_euler(), Lorene::Star_rot_Dirac::j_euler, Lorene::Scalar::set_etat_nondef(), Lorene::Tensor::set_etat_nondef(), and Lorene::Star_rot_Dirac::v2.
|
virtualinherited |
Ellipticity e.
Defined as
, where
and
are, respectively, the equatorial and polar circumferential radius, given by r_circ() and rp_circ().
Definition at line 431 of file strot_dirac_global.C.
References Lorene::Star_rot_Dirac::r_circ(), Lorene::Star_rot_Dirac::rp_circ(), and Lorene::sqrt().
|
inherited |
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition at line 465 of file star.C.
References Lorene::Param::add_int(), Lorene::Scalar::allocate_all(), Lorene::Star::del_deriv(), Lorene::Star::ener, Lorene::Eos::ener_ent(), Lorene::Star::ent, Lorene::Star::eos, Lorene::Mg3d::get_grille3d(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Star::mp, Lorene::Star::nbar, Lorene::Eos::nbar_ent(), Lorene::Star::nzet, Lorene::Star::press, Lorene::Eos::press_ent(), Lorene::Mtbl::set(), Lorene::Scalar::set_domain(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::std_spectral_base(), Lorene::Mtbl::t, and Lorene::Grille3d::x.
|
virtual |
Computes an equilibrium configuration.
| ent_c | [input] Central enthalpy |
| omega0 | [input] Requested central angular velocity (if fact_omega=1. ) |
| fact_omega | [input] 1.01 = search for the Keplerian frequency,
|
| nzadapt | [input] Number of (inner) domains where the mapping adaptation to an iso-enthalpy surface should be performed |
| ent_limit | [input] 1-D Tbl of dimension nzet which defines the enthalpy at the outer boundary of each domain |
| icontrol | [input] Set of integer parameters (stored as a 1-D Itbl of size 8) to control the iteration:
|
| control | [input] Set of parameters (stored as a 1-D Tbl of size 7) to control the iteration:
|
| mbar_wanted | [input] Requested baryon mass (effective only if mer_mass > mer_max ) |
| aexp_mass | [input] Exponent for the increase factor of the central enthalpy to converge to the requested baryon mass |
| diff | [output] 1-D Tbl of size 7 for the storage of some error indicators :
|
Definition at line 47 of file strot_dirac_diff_equil.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_t(), and Lorene::Star::mp.
|
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 101 of file star_equil_spher.C.
References Lorene::Star_rot::a_car, Lorene::Map_et::adapt(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Param::add_tbl(), Lorene::Scalar::annule(), Lorene::Star_rot::b_car, Lorene::Star_rot::bbb, Lorene::diffrel(), Lorene::Scalar::dsdr(), Lorene::Map_af::dsdr(), Lorene::Star_rot::dzeta, Lorene::Star::ener, Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::Star::equation_of_state(), Lorene::exp(), Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::Map_af::get_alpha(), Lorene::Map_et::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Map_et::get_beta(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map_af::homothetie(), Lorene::Scalar::integrale(), Lorene::Star::logn, Lorene::Star::mass_b(), Lorene::Star::mass_g(), Lorene::Star::mp, Lorene::Star::nn, Lorene::norme(), Lorene::Star::nzet, Lorene::Map_af::poisson(), Lorene::Star::press, Lorene::Star::s_euler, Lorene::Vector::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(), Lorene::Star::u_euler, Lorene::Scalar::val_grid_point(), and Lorene::Map::val_r().
| void Lorene::Star_rot_Dirac_diff::fait_omega_field | ( | double | omeg_min, |
| double | omeg_max, | ||
| double | precis, | ||
| int | nitermax | ||
| ) |
Computes
(member omega_field ).
The computation amounts to solving the equation
for
.
| omeg_min | [input] Lower bound of the interval for searching omega |
| omeg_max | [input] Higher bound of the interval for searching omega |
| precis | [input] Required precision in the determination of the zero by the secant method |
| nitermax | [input] Maximum number of iterations in the secant method to compute the zero. |
Definition at line 52 of file strot_dirac_diff_faitomeg.C.
References Lorene::Param::add_int(), Lorene::Param::add_scalar(), Lorene::Param::add_star(), Lorene::Scalar::annule(), Lorene::Star::beta, Lorene::Metric::cov(), Lorene::Scalar::div_rsint(), Lorene::Star::gamma, Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Star::mp, Lorene::Scalar::mult_rsint(), omega_field, Lorene::Star_rot_Dirac::psi4, and Lorene::Star_rot_Dirac::qqq.
| void Lorene::Star_rot_Dirac_diff::fait_prim_field | ( | ) |
Computes the member prim_field from omega_field .
Definition at line 176 of file strot_dirac_diff_faitomeg.C.
References Lorene::Scalar::allocate_all(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Star::mp, Lorene::Star::nzet, omega_field, par_frot, prim_field, primfrot, Lorene::Tbl::set(), Lorene::Scalar::set_domain(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().
| double Lorene::Star_rot_Dirac_diff::funct_omega | ( | double | omeg | ) | const |
Evaluates
, where F is the function defining the rotation profile.
This function is linked to the components of the fluid 4-velocity by
funct_omega calls frot with the parameters par_frot .
| omeg | [input] value of |
Definition at line 250 of file star_rot_dirac_diff.C.
|
inlineinherited |
|
inlineinherited |
Returns
.
Definition at line 255 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::aa_quad.
|
inlineinherited |
|
inlineinherited |
Returns the proper total energy density.
Definition at line 370 of file star.h.
References Lorene::Star::ener.
|
inlineinherited |
Returns the total energy density with respect to the Eulerian observer.
Definition at line 376 of file star.h.
References Lorene::Star::ener_euler.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
Returns the Lorentz factor between the fluid and Eulerian observers.
Definition at line 382 of file star.h.
References Lorene::Star::gam_euler.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
Returns the momentum density 3-vector with respect to the Eulerian observer.
Definition at line 227 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::j_euler.
|
inlineinherited |
Returns
.
Definition at line 217 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::ln_psi.
|
inlineinherited |
Returns the logarithm of the lapse N.
In the Newtonian case, this is the Newtonian gravitational potential (in units of
).
Definition at line 396 of file star.h.
References Lorene::Star::logn.
|
inlineinherited |
|
inlineinherited |
Returns the proper baryon density.
Definition at line 367 of file star.h.
References Lorene::Star::nbar.
|
inlineinherited |
|
inlineinherited |
Returns the number of domains occupied by the star.
Definition at line 358 of file star.h.
References Lorene::Star::nzet.
|
inlineinherited |
Returns the rotation angular velocity
.
Definition at line 196 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::omega.
|
virtual |
Returns the central value of the rotation angular velocity ([f_unit] )
Definition at line 262 of file star_rot_dirac_diff.C.
References omega_field, and Lorene::Scalar::val_grid_point().
|
inline |
Returns the angular velocity field
.
Definition at line 149 of file star_rot_dirac_diff.h.
References omega_field.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
Returns the conformal factor
.
Definition at line 202 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::psi4.
|
inlineinherited |
|
inlineinherited |
Returns the relativity parameter.
Definition at line 185 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::relat_type.
|
inlineinherited |
Returns the trace of the stress tensor in the Eulerian frame.
Definition at line 379 of file star.h.
References Lorene::Star::s_euler.
|
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 Lorene::Star::stress_euler.
|
inlineinherited |
|
inlineinherited |
Returns the conformal metric
.
Definition at line 240 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::tgamma.
|
inlineinherited |
Returns the fluid 3-velocity with respect to the Eulerian observer.
Definition at line 385 of file star.h.
References Lorene::Star::u_euler.
|
inlineinherited |
|
virtualinherited |
Error on the virial identity GRV2.
Definition at line 207 of file strot_dirac_global.C.
References Lorene::Star_rot_Dirac::aa, Lorene::Metric::con(), Lorene::Metric::cov(), Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::dsdr(), Lorene::Scalar::dsdt(), Lorene::Star::ener_euler, Lorene::Star::gamma, Lorene::Star::logn, Lorene::Star_rot_Dirac::p_grv2, Lorene::Star::press, Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), Lorene::Star_rot_Dirac::taa, and Lorene::Star_rot_Dirac::v2.
|
virtualinherited |
Error on the virial identity GRV3.
Definition at line 276 of file strot_dirac_global.C.
References Lorene::Star_rot_Dirac::aa_quad, Lorene::Metric::con(), Lorene::Metric::connect(), Lorene::contract(), Lorene::Scalar::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Metric::determinant(), Lorene::Star::gamma, Lorene::Connection::get_delta(), Lorene::Scalar::integrale(), Lorene::Star::logn, Lorene::Star_rot_Dirac::p_grv3, Lorene::Star::s_euler, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
virtual |
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame.
The calculation is performed starting from the quantities omega_field , ent , ener , press , and a_car ,
which are supposed to be up to date.
From these, the following fields are updated: gam_euler , u_euler , ener_euler , s_euler .
Reimplemented from Lorene::Star_rot_Dirac.
Definition at line 44 of file strot_dirac_diff_hydro.C.
References Lorene::Star::beta, Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Star_rot_Dirac::del_deriv(), Lorene::Star::ener, Lorene::Star::ener_euler, Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::Star_rot_Dirac::j_euler, Lorene::Scalar::mult_rsint(), Lorene::Star::nn, omega_field, Lorene::Star::press, Lorene::Star::s_euler, Lorene::Vector::set(), Lorene::Scalar::set_etat_zero(), Lorene::sqrt(), Lorene::Vector::std_spectral_base(), Lorene::Tensor::std_spectral_base(), Lorene::Scalar::std_spectral_base(), Lorene::Star::stress_euler, Lorene::Star::u_euler, and Lorene::Star_rot_Dirac::v2.
|
inlineinherited |
Checks whether the star is computed using a relativistic theory.
Definition at line 188 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::relat_type.
|
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 66 of file star_global.C.
References Lorene::Star::ent, Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::get_spectral_va(), Lorene::Star::mp, Lorene::Star::nzet, Lorene::Star::p_l_surf, and Lorene::Star::p_xi_surf.
|
virtualinherited |
Baryonic mass.
Implements Lorene::Star.
Definition at line 97 of file strot_dirac_global.C.
References Lorene::Metric::determinant(), 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().
|
virtualinherited |
Gravitational mass.
Implements Lorene::Star.
Definition at line 120 of file strot_dirac_global.C.
References Lorene::Star::beta, Lorene::contract(), Lorene::Metric::cov(), Lorene::Metric::determinant(), Lorene::Star::ener_euler, Lorene::Star::gamma, Lorene::Scalar::integrale(), Lorene::Star_rot_Dirac::j_euler, Lorene::Star::nn, Lorene::Star::p_mass_g, Lorene::Star::s_euler, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
| void Lorene::Star_rot_Dirac_diff::operator= | ( | const Star_rot_Dirac_diff & | star | ) |
Assignment to another Star_rot_Dirac_diff.
Definition at line 126 of file star_rot_dirac_diff.C.
References frot, omega_field, omega_max, omega_min, Lorene::Star_rot_Dirac::operator=(), par_frot, prim_field, and primfrot.
|
protectedvirtual |
Operator >> (virtual function called by the operator <<).
Reimplemented from Lorene::Star_rot_Dirac.
Definition at line 170 of file star_rot_dirac_diff.C.
References Lorene::Star_rot_Dirac::angu_mom(), Lorene::Star_rot_Dirac::aplat(), Lorene::Star_rot_Dirac::ellipt(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Star_rot_Dirac::grv2(), Lorene::Star_rot_Dirac::grv3(), Lorene::Star_rot_Dirac::mass_g(), Lorene::Star::mp, Lorene::Star::nzet, Lorene::Star_rot_Dirac::omega, omega_field, omega_max, omega_min, Lorene::Star::operator>>(), par_frot, Lorene::pow(), Lorene::Star_rot_Dirac::r_circ(), Lorene::Star::ray_eq(), Lorene::Star_rot_Dirac::rp_circ(), tsw(), and Lorene::Scalar::val_grid_point().
| double Lorene::Star_rot_Dirac_diff::prim_funct_omega | ( | double | omeg | ) | const |
Evaluates the primitive of
, where F is the function defining the rotation profile.
| omeg | [input] value of |
Definition at line 256 of file star_rot_dirac_diff.C.
|
virtualinherited |
Circumferential equatorial radius.
Definition at line 334 of file strot_dirac_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, Lorene::Star::nzet, and Lorene::Star_rot_Dirac::p_r_circ.
|
inherited |
Coordinate radius at
,
[r_unit].
Definition at line 111 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq.
|
inherited |
Coordinate radius at
,
[r_unit].
Definition at line 236 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq_3pis2.
|
inherited |
Coordinate radius at
,
[r_unit].
Definition at line 189 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq_pi.
|
inherited |
Coordinate radius at
,
[r_unit].
Definition at line 141 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq_pis2.
|
inherited |
Coordinate radius at
[r_unit].
Definition at line 281 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_pole.
|
virtualinherited |
Circumferential polar radius.
Definition at line 360 of file strot_dirac_global.C.
References Lorene::Scalar::annule(), Lorene::Scalar::annule_hard(), Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Metric::cov(), Lorene::Scalar::dsdt(), Lorene::Star::gamma, Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Star::l_surf(), Lorene::Star::mp, Lorene::Star::nzet, Lorene::Star_rot_Dirac::p_rp_circ, Lorene::Scalar::set_grid_point(), Lorene::Scalar::set_spectral_va(), Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), Lorene::Map::tet, Lorene::Valeur::val_point_jk(), Lorene::Map::val_r(), and Lorene::Star::xi_surf().
|
virtual |
Save in a file.
Reimplemented from Lorene::Star_rot_Dirac.
Definition at line 151 of file star_rot_dirac_diff.C.
References Lorene::fwrite_be(), omega_field, omega_max, omega_min, par_frot, Lorene::Star::sauve(), Lorene::Tbl::sauve(), and Lorene::Scalar::sauve().
|
protectedvirtualinherited |
Sets to 0x0 all the pointers on derived quantities.
Reimplemented from Lorene::Star.
Definition at line 251 of file star_rot_dirac.C.
References Lorene::Star_rot_Dirac::p_angu_mom, Lorene::Star_rot_Dirac::p_grv2, Lorene::Star_rot_Dirac::p_grv3, Lorene::Star_rot_Dirac::p_r_circ, Lorene::Star_rot_Dirac::p_rp_circ, and Lorene::Star_rot_Dirac::p_tsw.
|
inherited |
Assignment of the enthalpy field.
Definition at line 382 of file star.C.
References Lorene::Star::del_deriv(), Lorene::Star::ent, and Lorene::Star::equation_of_state().
|
inlineinherited |
|
inherited |
Solution of the tensor Poisson equation for rotating stars in Dirac gauge.
Definition at line 81 of file strot_dirac_solvehij.C.
References Lorene::Star_rot_Dirac::aa, Lorene::Star::beta, Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), 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::Scalar::derive_lie(), Lorene::Vector::divergence(), Lorene::Map::flat_met_spher(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Star_rot_Dirac::hh, Lorene::Tensor::inc_dzpuis(), Lorene::Scalar::inc_dzpuis(), Lorene::Star_rot_Dirac::ln_psi, Lorene::Star::mp, Lorene::Star::nn, Lorene::Vector::ope_killing_conf(), Lorene::Sym_tensor_trans::poisson(), Lorene::Star_rot_Dirac::psi2, Lorene::Star_rot_Dirac::psi4, Lorene::Star_rot_Dirac::qqq, Lorene::Star::s_euler, Lorene::Tensor::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Star::stress_euler, Lorene::Star_rot_Dirac::tgamma, and Lorene::Tensor::trace().
|
inherited |
Solution of the two scalar Poisson equations for rotating stars in Dirac gauge.
Definition at line 65 of file strot_dirac_solvenq.C.
References Lorene::Star::ener_euler, Lorene::Scalar::poisson(), Lorene::Star_rot_Dirac::psi4, and Lorene::Star::s_euler.
|
inherited |
Solution of the two scalar Poisson equations for rotating stars in Dirac gauge.
Definition at line 79 of file strot_dirac_solvenq.C.
References Lorene::Star_rot_Dirac::aa_quad, Lorene::contract(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Map::flat_met_spher(), Lorene::Map::get_bvect_spher(), Lorene::Star_rot_Dirac::hh, Lorene::Tensor::inc_dzpuis(), Lorene::Star_rot_Dirac::ln_psi, Lorene::Star::logn, Lorene::Star::mp, Lorene::Scalar::poisson(), Lorene::Star_rot_Dirac::psi4, Lorene::Star_rot_Dirac::tgamma, and Lorene::Tensor::up_down().
|
inherited |
Solution of the two scalar Poisson equations for rotating stars in Dirac gauge.
Definition at line 108 of file strot_dirac_solvenq.C.
References Lorene::Star_rot_Dirac::aa_quad, Lorene::contract(), Lorene::Metric::cov(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Map::flat_met_spher(), Lorene::Scalar::get_etat(), Lorene::Star_rot_Dirac::hh, Lorene::Tenseur::inc_dzpuis(), Lorene::Star_rot_Dirac::ln_psi, Lorene::Star::mp, Lorene::Star::nn, Lorene::Scalar::poisson(), Lorene::Star_rot_Dirac::psi2, Lorene::Star_rot_Dirac::psi4, Lorene::Star_rot_Dirac::qqq, Lorene::Star::s_euler, and Lorene::Star_rot_Dirac::tgamma.
|
inherited |
Solution of the shift equation for rotating stars in Dirac gauge.
Definition at line 57 of file strot_dirac_solveshift.C.
References Lorene::Star_rot_Dirac::aa, Lorene::Star::beta, Lorene::Metric::connect(), Lorene::contract(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Vector::div_free(), Lorene::Vector::divergence(), Lorene::Map::flat_met_spher(), Lorene::Connection::get_delta(), Lorene::Star_rot_Dirac::hh, Lorene::Tensor::inc_dzpuis(), Lorene::Star_rot_Dirac::j_euler, Lorene::Star_rot_Dirac::ln_psi, Lorene::Star::mp, Lorene::Star::nn, Lorene::Vector_divfree::poisson(), Lorene::Star_rot_Dirac::psi4, Lorene::Vector::set(), Lorene::Scalar::set_dzpuis(), and Lorene::Star_rot_Dirac::tgamma.
|
inlineinherited |
Returns the filtering order.
Definition at line 191 of file star_rot_dirac.h.
References Lorene::Star_rot_Dirac::spectral_filter.
|
virtual |
Ratio T/W.
Reimplemented from Lorene::Star_rot_Dirac.
Definition at line 48 of file strot_dirac_diff_global.C.
References Lorene::contract(), Lorene::Metric::cov(), Lorene::Metric::determinant(), Lorene::Star::ener, Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::Map::get_bvect_spher(), Lorene::Scalar::integrale(), Lorene::Star_rot_Dirac::j_euler, Lorene::Star_rot_Dirac::mass_g(), Lorene::Star::mp, omega_field, Lorene::Star_rot_Dirac::p_tsw, Lorene::Vector::set(), Lorene::Scalar::set_etat_zero(), Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
inherited |
Computes metric quantities from known potentials.
The calculation is performed starting from qqq, logn, shift, hh, which are supposed to be up to date. From these, the following fields are updated: nnn, psi4, psi2, ln_psi, tgamma, aa, taa, and aa_quad.
Definition at line 65 of file strot_dirac_upmetr.C.
References Lorene::Star_rot_Dirac::aa, Lorene::Star_rot_Dirac::aa_quad, Lorene::Star::beta, Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::contract(), Lorene::Star_rot_Dirac::del_deriv(), Lorene::Sym_tensor::derive_lie(), Lorene::exp(), Lorene::Star_rot_Dirac::flat, Lorene::Star::gamma, Lorene::Star_rot_Dirac::hh, Lorene::Tensor::inc_dzpuis(), Lorene::Star_rot_Dirac::ln_psi, Lorene::log(), Lorene::Star::logn, Lorene::Sym_tensor::longit_pot(), Lorene::Star::nn, Lorene::Vector::ope_killing_conf(), Lorene::Star_rot_Dirac::psi2, Lorene::Star_rot_Dirac::psi4, Lorene::Star_rot_Dirac::qqq, Lorene::Star_rot_Dirac::relat_type, Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), Lorene::Star_rot_Dirac::taa, Lorene::Star_rot_Dirac::tgamma, and Lorene::Tensor::up_down().
|
inherited |
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate
on the surface at the collocation points in
.
The stellar surface is defined as the location where the enthalpy (member ent) vanishes.
Definition at line 92 of file star_global.C.
References Lorene::Star::l_surf(), Lorene::Star::p_l_surf, and Lorene::Star::p_xi_surf.
|
friend |
|
protectedinherited |
Definition at line 97 of file star_rot_dirac.h.
|
protectedinherited |
Definition at line 99 of file star_rot_dirac.h.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
flat metric
(spherical components)
Definition at line 94 of file star_rot_dirac.h.
|
protected |
Function
defining the rotation profile.
This function is linked to the components of the fluid 4-velocity by
The first argument of frot must be
; the second argument contains the parameters, the first of which being necessarily the central value of
.
Definition at line 64 of file star_rot_dirac_diff.h.
|
protectedinherited |
|
protectedinherited |
is defined by
.
We impose the Dirac gauge
explicitly by defining
to be a symmetric transverse tensor.
Definition at line 106 of file star_rot_dirac.h.
|
protectedinherited |
Momentum density 3-vector with respect to the Eulerian observer.
Definition at line 87 of file star_rot_dirac.h.
|
protectedinherited |
Definition at line 77 of file star_rot_dirac.h.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Rotation angular velocity ([f_unit] )
Definition at line 69 of file star_rot_dirac.h.
|
protected |
Field
.
Definition at line 84 of file star_rot_dirac_diff.h.
|
protected |
Maximum value of
.
Definition at line 87 of file star_rot_dirac_diff.h.
|
protected |
Minimum value of
.
Definition at line 86 of file star_rot_dirac_diff.h.
|
mutableprotectedinherited |
Angular momentum.
Definition at line 116 of file star_rot_dirac.h.
|
mutableprotectedinherited |
Error on the virial identity GRV2.
Definition at line 117 of file star_rot_dirac.h.
|
mutableprotectedinherited |
Error on the virial identity GRV3.
Definition at line 118 of file star_rot_dirac.h.
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
Circumferential equatorial radius.
Definition at line 120 of file star_rot_dirac.h.
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
Circumferential polar radius.
Definition at line 121 of file star_rot_dirac.h.
|
mutableprotectedinherited |
Ratio T/W.
Definition at line 119 of file star_rot_dirac.h.
|
mutableprotectedinherited |
|
protected |
Parameters of the function
.
To be used as the second argument of functions frot and primfrot . The parameter par_frot(0) must always be the central angular velocity.
Definition at line 81 of file star_rot_dirac_diff.h.
|
protectedinherited |
|
protected |
Field
.
Definition at line 90 of file star_rot_dirac_diff.h.
|
protected |
Primitive of the function
, which vanishes at the stellar center.
The first argument of primfrot must be
; the second argument contains the parameters, the first of which being necessarily the central value of
.
Definition at line 72 of file star_rot_dirac_diff.h.
|
protectedinherited |
Definition at line 75 of file star_rot_dirac.h.
|
protectedinherited |
Conformal factor
.
Definition at line 74 of file star_rot_dirac.h.
|
protectedinherited |
Definition at line 76 of file star_rot_dirac.h.
|
protectedinherited |
Relativistic flag.
0 : Newtonian theory (not implemented) 1 : full General Relativity, maximal slicing and Dirac's Gauge 2 : Conformal Flatness Condition (CFC) 3 : CFC with
(see Cordero-Carrion et al. 2009).
Definition at line 59 of file star_rot_dirac.h.
|
protectedinherited |
|
protectedinherited |
Spectral exponential filtering order.
If 0, no filtering is done (see also Scalar::exponential_filter_r). Filtering is performed only in shells containing matter, i.e. for domain numbers l such that
nzet .
Definition at line 67 of file star_rot_dirac.h.
|
protectedinherited |
|
protectedinherited |
Definition at line 98 of file star_rot_dirac.h.
|
protectedinherited |
Definition at line 96 of file star_rot_dirac.h.
|
protectedinherited |
|
protectedinherited |
Definition at line 88 of file star_rot_dirac.h.