LORENE
|
Base class for stationary compact objects (under development). More...
#include <compobj.h>
Public Member Functions | |
Compobj (Map &map_i) | |
Standard constructor. More... | |
Compobj (const Compobj &) | |
Copy constructor. More... | |
Compobj (Map &map_i, FILE *) | |
Constructor from a file (see sauve(FILE* ) ). More... | |
virtual | ~Compobj () |
Destructor. More... | |
void | operator= (const Compobj &) |
Assignment to another Compobj. More... | |
Map & | set_mp () |
Read/write of the mapping. More... | |
const Map & | get_mp () const |
Returns the mapping. More... | |
const Scalar & | get_nn () const |
Returns the lapse function N . More... | |
const Vector & | get_beta () const |
Returns the shift vector . More... | |
const Metric & | get_gamma () const |
Returns the 3-metric . More... | |
const Scalar & | get_ener_euler () const |
Returns the total energy density E in the Eulerian frame. More... | |
const Vector & | get_mom_euler () const |
Returns the total 3-momentum density in the Eulerian frame. More... | |
const Sym_tensor & | get_stress_euler () const |
Returns the stress tensor with respect to the Eulerian observer. More... | |
const Sym_tensor & | get_kk () const |
Returns the extrinsic curvature tensor . More... | |
virtual void | sauve (FILE *) const |
Save in a file. More... | |
void | gyoto_data (const char *file_name) const |
Save in a file for GYOTO. More... | |
virtual void | extrinsic_curvature () |
Computation of the extrinsic curvature. More... | |
virtual double | adm_mass () const |
ADM mass (computed as a surface integral at spatial infinity) 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 ostream & | operator>> (ostream &) const |
Operator >> (virtual function called by the operator <<). More... | |
Protected Attributes | |
Map & | mp |
Mapping describing the coordinate system (r,theta,phi) More... | |
Scalar | nn |
Lapse function N . More... | |
Vector | beta |
Shift vector . More... | |
Metric | gamma |
3-metric More... | |
Scalar | ener_euler |
Total energy density E in the Eulerian frame. More... | |
Vector | mom_euler |
Total 3-momentum density in the Eulerian frame. More... | |
Sym_tensor | stress_euler |
Stress tensor with respect to the Eulerian observer. More... | |
Sym_tensor | kk |
Extrinsic curvature tensor . More... | |
double * | p_adm_mass |
ADM mass. More... | |
Friends | |
ostream & | operator<< (ostream &, const Compobj &) |
Display. More... | |
Base class for stationary compact objects (under development).
()
A Compobj
describes a single compact object (star or black hole), in a stationary state.
The spacetime metric is written according to the 3+1 formalism :
where is the 3-metric, described by a Lorene object of class Metric
.
The total energy-momentum tensor is orthogonally split with respect to the Eulerian observer as follows:
Lorene::Compobj::Compobj | ( | Map & | map_i | ) |
Standard constructor.
mp_i | Mapping on which the object is defined |
Definition at line 85 of file compobj.C.
References beta, ener_euler, kk, mom_euler, nn, set_der_0x0(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::std_spectral_base(), and stress_euler.
Lorene::Compobj::Compobj | ( | const Compobj & | co | ) |
Lorene::Compobj::Compobj | ( | Map & | map_i, |
FILE * | fich | ||
) |
Constructor from a file (see sauve(FILE* )
).
mp_i | Mapping on which the object is defined |
fich | input file (must have been created by the function sauve ) |
Definition at line 129 of file compobj.C.
References set_der_0x0().
|
virtual |
|
virtual |
ADM mass (computed as a surface integral at spatial infinity)
Definition at line 313 of file compobj.C.
References Lorene::Metric::cov(), Lorene::Tensor::derive_con(), Lorene::Tensor_sym::derive_con(), Lorene::Vector::flux(), gamma, Lorene::Tensor::get_triad(), mp, p_adm_mass, Lorene::Tensor::trace(), and Lorene::Tensor::up().
|
protectedvirtual |
Deletes all the derived quantities.
Reimplemented in Lorene::ScalarBH, Lorene::AltBH_QI, Lorene::Kerr_QI, Lorene::Star_QI, Lorene::Compobj_QI, and Lorene::Boson_star.
Definition at line 158 of file compobj.C.
References p_adm_mass, and set_der_0x0().
|
virtual |
Computation of the extrinsic curvature.
Reimplemented in Lorene::AltBH_QI, and Lorene::Compobj_QI.
Definition at line 293 of file compobj.C.
References beta, Lorene::Tensor::derive_cov(), Lorene::Tensor::down(), gamma, kk, nn, Lorene::Tensor::set(), and Lorene::Tensor::set_etat_qcq().
|
inline |
|
inline |
Returns the total energy density E in the Eulerian frame.
Definition at line 222 of file compobj.h.
References ener_euler.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Returns the stress tensor with respect to the Eulerian observer.
Definition at line 228 of file compobj.h.
References stress_euler.
void Lorene::Compobj::gyoto_data | ( | const char * | file_name | ) | const |
Save in a file for GYOTO.
Definition at line 213 of file compobj.C.
References beta, Lorene::Metric::con(), Lorene::Metric::cov(), Lorene::fwrite_be(), gamma, Lorene::Map::get_mg(), kk, mp, nn, Lorene::Mg3d::sauve(), Lorene::Map::sauve(), Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), and Lorene::Tensor_sym::sauve().
void Lorene::Compobj::operator= | ( | const Compobj & | co | ) |
Assignment to another Compobj.
Definition at line 178 of file compobj.C.
References beta, del_deriv(), ener_euler, gamma, kk, mom_euler, mp, nn, and stress_euler.
|
protectedvirtual |
Operator >> (virtual function called by the operator <<).
Reimplemented in Lorene::HiggsMonopole, Lorene::ScalarBH, Lorene::AltBH_QI, Lorene::Kerr_QI, Lorene::Star_QI, Lorene::Compobj_QI, and Lorene::Boson_star.
Definition at line 242 of file compobj.C.
References Lorene::Metric::cov(), ener_euler, gamma, kk, mp, nn, stress_euler, and Lorene::Scalar::val_grid_point().
|
virtual |
Save in a file.
Reimplemented in Lorene::ScalarBH, Lorene::AltBH_QI, Lorene::Kerr_QI, Lorene::Star_QI, Lorene::Compobj_QI, and Lorene::Boson_star.
Definition at line 199 of file compobj.C.
References beta, ener_euler, gamma, mom_euler, nn, Lorene::Metric::sauve(), Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), Lorene::Tensor_sym::sauve(), and stress_euler.
|
protected |
Sets to 0x0
all the pointers on derived quantities.
Definition at line 166 of file compobj.C.
References p_adm_mass.
|
inline |
|
friend |
|
protected |
|
protected |
|
protected |
|
protected |
|
mutableprotected |
|
protected |