LORENE
Lorene::Star_rot_CFC Class Reference

Class for relativistic rotating stars in Conformal Flatness Condition and maximal slicing. More...

#include <star_rot_cfc.h>

Inheritance diagram for Lorene::Star_rot_CFC:
Lorene::Star

Public Member Functions

 Star_rot_CFC (Map &mp_i, int nzet_i, const Eos &eos_i, int relat_i=3, int filter=0)
 Standard constructor. More...
 
 Star_rot_CFC (const Star_rot_CFC &)
 Copy constructor. More...
 
 Star_rot_CFC (Map &mp_i, const Eos &eos_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ). More...
 
virtual ~Star_rot_CFC ()
 Destructor. More...
 
void operator= (const Star_rot_CFC &)
 Assignment to another Star_rot_CFC. 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 $\Omega$. More...
 
const Scalarget_psi4 () const
 Returns the conformal factor $\Psi^4$. More...
 
const Scalarget_psi2 () const
 Returns $\Psi^2$. More...
 
const Scalarget_psi () const
 Returns $\Psi$. More...
 
const Vectorget_j_euler () const
 Returns the momentum density 3-vector with respect to the Eulerian observer. More...
 
const Scalarget_v2 () const
 Returns $\gamma_{ij}v^i v^j$. More...
 
const Sym_tensor get_hatA () const
 Returns $\hat{A}^{ij}$. More...
 
const Scalar get_hatA_quad () const
 Returns $\tilde{A}_{ij} A^{ij}$. More...
 
virtual void sauve (FILE *) const
 Save in a file. 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 tsw () const
 Ratio T/W. 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...
 
virtual void hydro_euler ()
 Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame. More...
 
void update_metric ()
 Computes metric quantities from known potentials. More...
 
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...
 
void solve_logn_f (Scalar &ln_f_new) const
 Solution of the `‘matter’' part of the Poisson equation for the lapse for rotating stars in CFC. More...
 
void solve_logn_q (Scalar &ln_q_new) const
 Solution of the quadratic part of the Poisson equation for the lapse for rotating stars in CFC. More...
 
void solve_psi (Scalar &psi_new)
 Solution of the equations for the conformal factor for rotating stars in CFC. More...
 
void solve_shift (Vector &shift_new) const
 Solution of the shift equation for rotating stars in CFC. More...
 
Mapset_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 Mapget_mp () const
 Returns the mapping. More...
 
int get_nzet () const
 Returns the number of domains occupied by the star. More...
 
const Eosget_eos () const
 Returns the equation of state. More...
 
const Scalarget_ent () const
 Returns the enthalpy field. More...
 
const Scalarget_nbar () const
 Returns the proper baryon density. More...
 
const Scalarget_ener () const
 Returns the proper total energy density. More...
 
const Scalarget_press () const
 Returns the fluid pressure. More...
 
const Scalarget_ener_euler () const
 Returns the total energy density with respect to the Eulerian observer. More...
 
const Scalarget_s_euler () const
 Returns the trace of the stress tensor in the Eulerian frame. More...
 
const Scalarget_gam_euler () const
 Returns the Lorentz factor between the fluid and Eulerian observers. More...
 
const Vectorget_u_euler () const
 Returns the fluid 3-velocity with respect to the Eulerian observer. More...
 
const Tensorget_stress_euler () const
 Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer. More...
 
const Scalarget_logn () const
 Returns the logarithm of the lapse N. More...
 
const Scalarget_nn () const
 Returns the lapse function N. More...
 
const Vectorget_beta () const
 Returns the shift vector $\beta^i$. More...
 
const Scalarget_lnq () const
 
const Metricget_gamma () const
 Returns the 3-metric $\gamma$. More...
 
double ray_eq () const
 Coordinate radius at $\phi=0$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq_pis2 () const
 Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq_pi () const
 Coordinate radius at $\phi=\pi$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq_3pis2 () const
 Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_pole () const
 Coordinate radius at $\theta=0$ [r_unit]. More...
 
virtual const Itbll_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 $(\theta', \phi')$. More...
 
const Tblxi_surf () const
 Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$. More...
 

Protected Member Functions

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...
 
virtual ostream & operator>> (ostream &) const
 Operator >> (virtual function called by the operator <<). More...
 

Protected Attributes

int relat_type
 Relativistic flag. More...
 
int spectral_filter
 Spectral exponential filtering order. More...
 
double omega
 Rotation angular velocity ([f_unit] ) More...
 
Scalar psi
 Conformal factor $\Psi$. More...
 
Vector j_euler
 Momentum density 3-vector with respect to the Eulerian observer. More...
 
Scalar v2
 $\gamma_{ij}v^i v^j$ More...
 
Scalar psi4
 Conformal factor $\Psi^4$. More...
 
Scalar psi2
 $\Psi^2$ More...
 
const Metric_flatflat
 flat metric $f_{ij}$ (spherical components) More...
 
Sym_tensor hatA
 $\hat{A}^{ij}$ More...
 
Scalar hatA_quad
 $f_{il}\, f_{jm}\, \hat{A}_{ij} \hat{A}^{lm}$ 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...
 
Mapmp
 Mapping associated with the star. More...
 
int nzet
 Number of domains of *mp occupied by the star. More...
 
const Eoseos
 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 $\phi=0$, $\theta=\pi/2$. More...
 
double * p_ray_eq_pis2
 Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$. More...
 
double * p_ray_eq_pi
 Coordinate radius at $\phi=\pi$, $\theta=\pi/2$. More...
 
double * p_ray_eq_3pis2
 Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$. More...
 
double * p_ray_pole
 Coordinate radius at $\theta=0$. More...
 
Itblp_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 $(\theta', \phi')$. More...
 
Tblp_xi_surf
 Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$. More...
 
double * p_mass_b
 Baryon mass. More...
 
double * p_mass_g
 Gravitational mass. More...
 

Detailed Description

Class for relativistic rotating stars in Conformal Flatness Condition and maximal slicing.

()

Definition at line 39 of file star_rot_cfc.h.

Constructor & Destructor Documentation

◆ Star_rot_CFC() [1/3]

Lorene::Star_rot_CFC::Star_rot_CFC ( Map mp_i,
int  nzet_i,
const Eos eos_i,
int  relat_i = 3,
int  filter = 0 
)

Standard constructor.

Parameters
mp_iMapping on which the star will be defined
nzet_iNumber of domains occupied by the star
eos_iEquation of state of the stellar matter
relat_iReltivity parameter (see relat_type )
filterorder for spectral exponential filtering

Definition at line 44 of file star_rot_cfc.C.

References hatA, hatA_quad, j_euler, omega, psi, psi2, psi4, relat_type, set_der_0x0(), Lorene::Tensor::set_etat_zero(), spectral_filter, Lorene::Scalar::std_spectral_base(), and v2.

◆ Star_rot_CFC() [2/3]

Lorene::Star_rot_CFC::Star_rot_CFC ( const Star_rot_CFC star)

Copy constructor.

Definition at line 83 of file star_rot_cfc.C.

References omega, and set_der_0x0().

◆ Star_rot_CFC() [3/3]

Lorene::Star_rot_CFC::Star_rot_CFC ( Map mp_i,
const Eos eos_i,
FILE *  fich 
)

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

Parameters
mp_iMapping on which the star will be defined
eos_iEquation of state of the stellar matter
fichinput file (must have been created by the function sauve )

Definition at line 107 of file star_rot_cfc.C.

References Lorene::Star::beta, Lorene::contract(), Lorene::Star::equation_of_state(), flat, Lorene::fread_be(), Lorene::Map::get_bvect_spher(), hatA, hatA_quad, hydro_euler(), j_euler, omega, Lorene::Vector::ope_killing_conf(), Lorene::Vector::poisson(), psi2, psi4, relat_type, set_der_0x0(), spectral_filter, Lorene::Tensor::up_down(), and update_metric().

◆ ~Star_rot_CFC()

Lorene::Star_rot_CFC::~Star_rot_CFC ( )
virtual

Destructor.

Definition at line 153 of file star_rot_cfc.C.

References del_deriv().

Member Function Documentation

◆ angu_mom()

◆ aplat()

double Lorene::Star_rot_CFC::aplat ( ) const
virtual

Flattening r_pole/r_eq.

Definition at line 369 of file strot_cfc_global.C.

References Lorene::Star::ray_eq(), and Lorene::Star::ray_pole().

◆ del_deriv()

void Lorene::Star_rot_CFC::del_deriv ( ) const
protectedvirtual

Deletes all the derived quantities.

Reimplemented from Lorene::Star.

Definition at line 164 of file star_rot_cfc.C.

References Lorene::Star::del_deriv(), p_angu_mom, p_grv2, p_grv3, p_r_circ, p_rp_circ, p_tsw, and set_der_0x0().

◆ del_hydro_euler()

void Lorene::Star_rot_CFC::del_hydro_euler ( )
protectedvirtual

Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.

Reimplemented from Lorene::Star.

Definition at line 192 of file star_rot_cfc.C.

References del_deriv(), Lorene::Star::del_hydro_euler(), j_euler, Lorene::Scalar::set_etat_nondef(), Lorene::Tensor::set_etat_nondef(), and v2.

◆ ellipt()

double Lorene::Star_rot_CFC::ellipt ( ) const
virtual

Ellipticity e.

Defined as $ e = \sqrt{1 - \left( \frac{R^c_e}{R^c_p} \right)^2} $, where $R^c_e$ and $R^c_p$ are, respectively, the equatorial and polar circumferential radius, given by r_circ() and rp_circ().

Definition at line 379 of file strot_cfc_global.C.

References r_circ(), rp_circ(), and Lorene::sqrt().

◆ equation_of_state()

◆ equilibrium()

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

Definition at line 36 of file strot_cfc_equil.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_t(), and Lorene::Star::mp.

◆ equilibrium_spher()

void Lorene::Star::equilibrium_spher ( double  ent_c,
double  precis = 1.e-14,
const Tbl pent_limit = 0x0 
)
virtualinherited

Computes a spherical static configuration.

Parameters
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().

◆ get_beta()

const Vector& Lorene::Star::get_beta ( ) const
inlineinherited

Returns the shift vector $\beta^i$.

Definition at line 402 of file star.h.

References Lorene::Star::beta.

◆ get_ener()

const Scalar& Lorene::Star::get_ener ( ) const
inlineinherited

Returns the proper total energy density.

Definition at line 370 of file star.h.

References Lorene::Star::ener.

◆ get_ener_euler()

const Scalar& Lorene::Star::get_ener_euler ( ) const
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.

◆ get_ent()

const Scalar& Lorene::Star::get_ent ( ) const
inlineinherited

Returns the enthalpy field.

Definition at line 364 of file star.h.

References Lorene::Star::ent.

◆ get_eos()

const Eos& Lorene::Star::get_eos ( ) const
inlineinherited

Returns the equation of state.

Definition at line 361 of file star.h.

References Lorene::Star::eos.

◆ get_gam_euler()

const Scalar& Lorene::Star::get_gam_euler ( ) const
inlineinherited

Returns the Lorentz factor between the fluid and Eulerian observers.

Definition at line 382 of file star.h.

References Lorene::Star::gam_euler.

◆ get_gamma()

const Metric& Lorene::Star::get_gamma ( ) const
inlineinherited

Returns the 3-metric $\gamma$.

Definition at line 409 of file star.h.

References Lorene::Star::gamma.

◆ get_hatA()

const Sym_tensor Lorene::Star_rot_CFC::get_hatA ( ) const
inline

Returns $\hat{A}^{ij}$.

Definition at line 215 of file star_rot_cfc.h.

References hatA.

◆ get_hatA_quad()

const Scalar Lorene::Star_rot_CFC::get_hatA_quad ( ) const
inline

Returns $\tilde{A}_{ij} A^{ij}$.

Definition at line 220 of file star_rot_cfc.h.

References hatA_quad.

◆ get_j_euler()

const Vector& Lorene::Star_rot_CFC::get_j_euler ( ) const
inline

Returns the momentum density 3-vector with respect to the Eulerian observer.

Definition at line 202 of file star_rot_cfc.h.

References j_euler.

◆ get_logn()

const Scalar& Lorene::Star::get_logn ( ) const
inlineinherited

Returns the logarithm of the lapse N.

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 396 of file star.h.

References Lorene::Star::logn.

◆ get_mp()

const Map& Lorene::Star::get_mp ( ) const
inlineinherited

Returns the mapping.

Definition at line 355 of file star.h.

References Lorene::Star::mp.

◆ get_nbar()

const Scalar& Lorene::Star::get_nbar ( ) const
inlineinherited

Returns the proper baryon density.

Definition at line 367 of file star.h.

References Lorene::Star::nbar.

◆ get_nn()

const Scalar& Lorene::Star::get_nn ( ) const
inlineinherited

Returns the lapse function N.

Definition at line 399 of file star.h.

References Lorene::Star::nn.

◆ get_nzet()

int Lorene::Star::get_nzet ( ) const
inlineinherited

Returns the number of domains occupied by the star.

Definition at line 358 of file star.h.

References Lorene::Star::nzet.

◆ get_omega()

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

Returns the rotation angular velocity $\Omega$.

Definition at line 177 of file star_rot_cfc.h.

References omega.

◆ get_press()

const Scalar& Lorene::Star::get_press ( ) const
inlineinherited

Returns the fluid pressure.

Definition at line 373 of file star.h.

References Lorene::Star::press.

◆ get_psi()

const Scalar& Lorene::Star_rot_CFC::get_psi ( ) const
inline

Returns $\Psi$.

Definition at line 193 of file star_rot_cfc.h.

References psi.

◆ get_psi2()

const Scalar& Lorene::Star_rot_CFC::get_psi2 ( ) const
inline

Returns $\Psi^2$.

Definition at line 188 of file star_rot_cfc.h.

References psi2.

◆ get_psi4()

const Scalar& Lorene::Star_rot_CFC::get_psi4 ( ) const
inline

Returns the conformal factor $\Psi^4$.

Definition at line 183 of file star_rot_cfc.h.

References psi4.

◆ get_relat()

int Lorene::Star_rot_CFC::get_relat ( ) const
inline

Returns the relativity parameter.

Definition at line 166 of file star_rot_cfc.h.

References relat_type.

◆ get_s_euler()

const Scalar& Lorene::Star::get_s_euler ( ) const
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.

◆ get_stress_euler()

const Tensor& Lorene::Star::get_stress_euler ( ) const
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.

◆ get_u_euler()

const Vector& Lorene::Star::get_u_euler ( ) const
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.

◆ get_v2()

const Scalar& Lorene::Star_rot_CFC::get_v2 ( ) const
inline

Returns $\gamma_{ij}v^i v^j$.

Definition at line 207 of file star_rot_cfc.h.

References v2.

◆ grv2()

◆ grv3()

◆ hydro_euler()

◆ is_relativistic()

bool Lorene::Star_rot_CFC::is_relativistic ( ) const
inline

Checks whether the star is computed using a relativistic theory.

Definition at line 169 of file star_rot_cfc.h.

References relat_type.

◆ l_surf()

const Itbl & Lorene::Star::l_surf ( ) const
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 $(\theta', \phi')$.

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.

◆ mass_b()

double Lorene::Star_rot_CFC::mass_b ( ) const
virtual

◆ mass_g()

◆ operator=()

void Lorene::Star_rot_CFC::operator= ( const Star_rot_CFC star)

Assignment to another Star_rot_CFC.

Definition at line 212 of file star_rot_cfc.C.

References del_deriv(), flat, hatA, hatA_quad, j_euler, omega, Lorene::Star::operator=(), psi, psi2, psi4, relat_type, spectral_filter, and v2.

◆ operator>>()

ostream & Lorene::Star_rot_CFC::operator>> ( ostream &  ost) const
protectedvirtual

◆ r_circ()

double Lorene::Star_rot_CFC::r_circ ( ) const
virtual

Circumferential equatorial radius.

Definition at line 282 of file strot_cfc_global.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, Lorene::Star::nzet, and p_r_circ.

◆ ray_eq()

double Lorene::Star::ray_eq ( ) const
inherited

◆ ray_eq_3pis2()

double Lorene::Star::ray_eq_3pis2 ( ) const
inherited

◆ ray_eq_pi()

double Lorene::Star::ray_eq_pi ( ) const
inherited

◆ ray_eq_pis2()

double Lorene::Star::ray_eq_pis2 ( ) const
inherited

◆ ray_pole()

double Lorene::Star::ray_pole ( ) const
inherited

Coordinate radius at $\theta=0$ [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.

◆ rp_circ()

◆ sauve()

void Lorene::Star_rot_CFC::sauve ( FILE *  fich) const
virtual

◆ set_der_0x0()

void Lorene::Star_rot_CFC::set_der_0x0 ( ) const
protectedvirtual

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Lorene::Star.

Definition at line 180 of file star_rot_cfc.C.

References p_angu_mom, p_grv2, p_grv3, p_r_circ, p_rp_circ, and p_tsw.

◆ set_enthalpy()

void Lorene::Star::set_enthalpy ( const Scalar ent_i)
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().

◆ set_mp()

Map& Lorene::Star::set_mp ( )
inlineinherited

Read/write of the mapping.

Definition at line 322 of file star.h.

References Lorene::Star::mp.

◆ solve_logn_f()

void Lorene::Star_rot_CFC::solve_logn_f ( Scalar ln_f_new) const

Solution of the `‘matter’' part of the Poisson equation for the lapse for rotating stars in CFC.

Definition at line 39 of file strot_cfc_metric.C.

References Lorene::Star::ener_euler, Lorene::Scalar::poisson(), psi4, and Lorene::Star::s_euler.

◆ solve_logn_q()

void Lorene::Star_rot_CFC::solve_logn_q ( Scalar ln_q_new) const

Solution of the quadratic part of the Poisson equation for the lapse for rotating stars in CFC.

Definition at line 53 of file strot_cfc_metric.C.

References Lorene::contract(), Lorene::Scalar::derive_con(), Lorene::Scalar::derive_cov(), flat, hatA_quad, Lorene::log(), Lorene::Star::logn, Lorene::Scalar::poisson(), psi, psi4, Lorene::Scalar::std_spectral_base(), and Lorene::Tensor::up_down().

◆ solve_psi()

void Lorene::Star_rot_CFC::solve_psi ( Scalar psi_new)

◆ solve_shift()

void Lorene::Star_rot_CFC::solve_shift ( Vector shift_new) const

Solution of the shift equation for rotating stars in CFC.

Definition at line 106 of file strot_cfc_metric.C.

References Lorene::Sym_tensor::divergence(), flat, hatA, Lorene::Tensor::inc_dzpuis(), Lorene::Star::nn, Lorene::Vector::poisson(), psi2, psi4, and Lorene::Vector::set().

◆ spectral_filter_order()

int Lorene::Star_rot_CFC::spectral_filter_order ( ) const
inline

Returns the filtering order.

Definition at line 172 of file star_rot_cfc.h.

References spectral_filter.

◆ tsw()

double Lorene::Star_rot_CFC::tsw ( ) const
virtual

◆ update_metric()

void Lorene::Star_rot_CFC::update_metric ( )

Computes metric quantities from known potentials.

The calculation is performed starting from psi, logn, hatA, which are supposed to be up to date. From these, the following fields are updated: nnn, psi4, psi2,shift, and hatA_quad.

Definition at line 125 of file strot_cfc_metric.C.

References Lorene::Metric_flat::cov(), del_deriv(), Lorene::exp(), flat, Lorene::Star::gamma, Lorene::Star::logn, Lorene::Star::nn, psi, psi2, psi4, and Lorene::Scalar::std_spectral_base().

◆ xi_surf()

const Tbl & Lorene::Star::xi_surf ( ) const
inherited

Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.

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.

Member Data Documentation

◆ beta

Vector Lorene::Star::beta
protectedinherited

Shift vector.

Definition at line 228 of file star.h.

◆ ener

Scalar Lorene::Star::ener
protectedinherited

Total energy density in the fluid frame.

Definition at line 193 of file star.h.

◆ ener_euler

Scalar Lorene::Star::ener_euler
protectedinherited

Total energy density in the Eulerian frame.

Definition at line 198 of file star.h.

◆ ent

Scalar Lorene::Star::ent
protectedinherited

Log-enthalpy.

Definition at line 190 of file star.h.

◆ eos

const Eos& Lorene::Star::eos
protectedinherited

Equation of state of the stellar matter.

Definition at line 185 of file star.h.

◆ flat

const Metric_flat& Lorene::Star_rot_CFC::flat
protected

flat metric $f_{ij}$ (spherical components)

Definition at line 85 of file star_rot_cfc.h.

◆ gam_euler

Scalar Lorene::Star::gam_euler
protectedinherited

Lorentz factor between the fluid and Eulerian observers.

Definition at line 204 of file star.h.

◆ gamma

Metric Lorene::Star::gamma
protectedinherited

3-metric

Definition at line 235 of file star.h.

◆ hatA

Sym_tensor Lorene::Star_rot_CFC::hatA
protected

$\hat{A}^{ij}$

Definition at line 87 of file star_rot_cfc.h.

◆ hatA_quad

Scalar Lorene::Star_rot_CFC::hatA_quad
protected

$f_{il}\, f_{jm}\, \hat{A}_{ij} \hat{A}^{lm}$

Definition at line 88 of file star_rot_cfc.h.

◆ j_euler

Vector Lorene::Star_rot_CFC::j_euler
protected

Momentum density 3-vector with respect to the Eulerian observer.

Definition at line 75 of file star_rot_cfc.h.

◆ logn

Scalar Lorene::Star::logn
protectedinherited

Logarithm of the lapse N .

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 222 of file star.h.

◆ mp

Map& Lorene::Star::mp
protectedinherited

Mapping associated with the star.

Definition at line 180 of file star.h.

◆ nbar

Scalar Lorene::Star::nbar
protectedinherited

Baryon density in the fluid frame.

Definition at line 192 of file star.h.

◆ nn

Scalar Lorene::Star::nn
protectedinherited

Lapse function N .

Definition at line 225 of file star.h.

◆ nzet

int Lorene::Star::nzet
protectedinherited

Number of domains of *mp occupied by the star.

Definition at line 183 of file star.h.

◆ omega

double Lorene::Star_rot_CFC::omega
protected

Rotation angular velocity ([f_unit] )

Definition at line 60 of file star_rot_cfc.h.

◆ p_angu_mom

double* Lorene::Star_rot_CFC::p_angu_mom
mutableprotected

Angular momentum.

Definition at line 97 of file star_rot_cfc.h.

◆ p_grv2

double* Lorene::Star_rot_CFC::p_grv2
mutableprotected

Error on the virial identity GRV2.

Definition at line 98 of file star_rot_cfc.h.

◆ p_grv3

double* Lorene::Star_rot_CFC::p_grv3
mutableprotected

Error on the virial identity GRV3.

Definition at line 99 of file star_rot_cfc.h.

◆ p_l_surf

Itbl* Lorene::Star::p_l_surf
mutableprotectedinherited

Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$.

Definition at line 260 of file star.h.

◆ p_mass_b

double* Lorene::Star::p_mass_b
mutableprotectedinherited

Baryon mass.

Definition at line 268 of file star.h.

◆ p_mass_g

double* Lorene::Star::p_mass_g
mutableprotectedinherited

Gravitational mass.

Definition at line 269 of file star.h.

◆ p_r_circ

double* Lorene::Star_rot_CFC::p_r_circ
mutableprotected

Circumferential equatorial radius.

Definition at line 101 of file star_rot_cfc.h.

◆ p_ray_eq

double* Lorene::Star::p_ray_eq
mutableprotectedinherited

Coordinate radius at $\phi=0$, $\theta=\pi/2$.

Definition at line 242 of file star.h.

◆ p_ray_eq_3pis2

double* Lorene::Star::p_ray_eq_3pis2
mutableprotectedinherited

Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$.

Definition at line 251 of file star.h.

◆ p_ray_eq_pi

double* Lorene::Star::p_ray_eq_pi
mutableprotectedinherited

Coordinate radius at $\phi=\pi$, $\theta=\pi/2$.

Definition at line 248 of file star.h.

◆ p_ray_eq_pis2

double* Lorene::Star::p_ray_eq_pis2
mutableprotectedinherited

Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$.

Definition at line 245 of file star.h.

◆ p_ray_pole

double* Lorene::Star::p_ray_pole
mutableprotectedinherited

Coordinate radius at $\theta=0$.

Definition at line 254 of file star.h.

◆ p_rp_circ

double* Lorene::Star_rot_CFC::p_rp_circ
mutableprotected

Circumferential polar radius.

Definition at line 102 of file star_rot_cfc.h.

◆ p_tsw

double* Lorene::Star_rot_CFC::p_tsw
mutableprotected

Ratio T/W.

Definition at line 100 of file star_rot_cfc.h.

◆ p_xi_surf

Tbl* Lorene::Star::p_xi_surf
mutableprotectedinherited

Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.

Definition at line 266 of file star.h.

◆ press

Scalar Lorene::Star::press
protectedinherited

Fluid pressure.

Definition at line 194 of file star.h.

◆ psi

Scalar Lorene::Star_rot_CFC::psi
protected

Conformal factor $\Psi$.

Definition at line 65 of file star_rot_cfc.h.

◆ psi2

Scalar Lorene::Star_rot_CFC::psi2
protected

$\Psi^2$

Definition at line 82 of file star_rot_cfc.h.

◆ psi4

Scalar Lorene::Star_rot_CFC::psi4
protected

Conformal factor $\Psi^4$.

Definition at line 81 of file star_rot_cfc.h.

◆ relat_type

int Lorene::Star_rot_CFC::relat_type
protected

Relativistic flag.

0 : Newtonian theory (not implemented) 3 : CFC with $\hat{A}^{ij}_{TT} = 0$ (see Cordero-Carrion et al. 2009).

Definition at line 50 of file star_rot_cfc.h.

◆ s_euler

Scalar Lorene::Star::s_euler
protectedinherited

Trace of the stress scalar in the Eulerian frame.

Definition at line 201 of file star.h.

◆ spectral_filter

int Lorene::Star_rot_CFC::spectral_filter
protected

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 $0 < l <$ nzet .

Definition at line 58 of file star_rot_cfc.h.

◆ stress_euler

Sym_tensor Lorene::Star::stress_euler
protectedinherited

Spatial part of the stress-energy tensor with respect to the Eulerian observer.

Definition at line 212 of file star.h.

◆ u_euler

Vector Lorene::Star::u_euler
protectedinherited

Fluid 3-velocity with respect to the Eulerian observer.

Definition at line 207 of file star.h.

◆ v2

Scalar Lorene::Star_rot_CFC::v2
protected

$\gamma_{ij}v^i v^j$

Definition at line 76 of file star_rot_cfc.h.


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