LORENE

Base class for stationary compact objects (under development). More...

#include <compobj.h>

Inheritance diagram for Lorene::Compobj:
Lorene::Compobj_QI Lorene::HiggsMonopole Lorene::ScalarBH Lorene::AltBH_QI Lorene::Kerr_QI Lorene::Star_QI Lorene::Boson_star

Public Member Functions

 Compobj (Map &map_i)
 Standard constructor.
 Compobj (const Compobj &)
 Copy constructor.
 Compobj (Map &map_i, FILE *)
 Constructor from a file (see sauve(FILE* )).
virtual ~Compobj ()
 Destructor.
void operator= (const Compobj &)
 Assignment to another Compobj.
Mapset_mp ()
 Read/write of the mapping.
const Mapget_mp () const
 Returns the mapping.
const Scalarget_nn () const
 Returns the lapse function N .
const Vectorget_beta () const
 Returns the shift vector $\beta^i$.
const Metricget_gamma () const
 Returns the 3-metric $\gamma_{ij}$.
const Scalarget_ener_euler () const
 Returns the total energy density E in the Eulerian frame.
const Vectorget_mom_euler () const
 Returns the total 3-momentum density $P^i$ in the Eulerian frame.
const Sym_tensorget_stress_euler () const
 Returns the stress tensor $S_{ij}$ with respect to the Eulerian observer.
const Sym_tensorget_kk () const
 Returns the extrinsic curvature tensor $K_{ij}$.
virtual void sauve (FILE *) const
 Save in a file.
void gyoto_data (const char *file_name) const
 Save in a file for GYOTO.
virtual void extrinsic_curvature ()
 Computation of the extrinsic curvature.
virtual double adm_mass () const
 ADM mass (computed as a surface integral at spatial infinity).

Protected Member Functions

virtual void del_deriv () const
 Deletes all the derived quantities.
void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities.
virtual ostream & operator>> (ostream &) const
 Operator >> (virtual function called by the operator <<).

Protected Attributes

Mapmp
 Mapping describing the coordinate system (r,theta,phi).
Scalar nn
 Lapse function N .
Vector beta
 Shift vector $\beta^i$.
Metric gamma
 3-metric $\gamma_{ij}$
Scalar ener_euler
 Total energy density E in the Eulerian frame.
Vector mom_euler
 Total 3-momentum density $P^i$ in the Eulerian frame.
Sym_tensor stress_euler
 Stress tensor $S_{ij}$ with respect to the Eulerian observer.
Sym_tensor kk
 Extrinsic curvature tensor $K_{ij}$.
double * p_adm_mass
 ADM mass.

Friends

ostream & operator<< (ostream &ost, const Compobj &co)
 Display.

Detailed Description

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 :

\‍[  ds^2 = - N^2  dt^2 + \gamma_{ij} ( dx^i + \beta^i dt )
              (dx^j + \beta^j dt )
\‍]

where $\gamma_{ij}$ 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:

\‍[   T_{\alpha\beta} = E n_\alpha n_\beta + P_\alpha n_\beta + n_\alpha P_\beta + S_{\alpha\beta}
\‍]

Definition at line 126 of file compobj.h.

Constructor & Destructor Documentation

◆ Compobj() [1/3]

Lorene::Compobj::Compobj ( Map & map_i)

Standard constructor.

Parameters
mp_iMapping on which the object is defined

Definition at line 82 of file compobj.C.

References beta, ener_euler, Lorene::flat_met_spher(), gamma, Lorene::get_bvect_spher(), kk, Lorene::Map(), mom_euler, mp, nn, set_der_0x0(), and stress_euler.

◆ Compobj() [2/3]

Lorene::Compobj::Compobj ( const Compobj & co)

Copy constructor.

Definition at line 109 of file compobj.C.

References beta, Compobj(), ener_euler, gamma, kk, mom_euler, mp, nn, set_der_0x0(), and stress_euler.

◆ Compobj() [3/3]

Lorene::Compobj::Compobj ( Map & map_i,
FILE * fich )

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

Parameters
mp_iMapping on which the object is defined
fichinput file (must have been created by the function sauve)

Definition at line 126 of file compobj.C.

References beta, ener_euler, gamma, Lorene::get_bvect_spher(), Lorene::get_mg(), kk, Lorene::Map(), mom_euler, mp, nn, set_der_0x0(), and stress_euler.

◆ ~Compobj()

Lorene::Compobj::~Compobj ( )
virtual

Destructor.

Definition at line 144 of file compobj.C.

References del_deriv().

Member Function Documentation

◆ adm_mass()

double Lorene::Compobj::adm_mass ( ) const
virtual

ADM mass (computed as a surface integral at spatial infinity).

Definition at line 310 of file compobj.C.

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

◆ del_deriv()

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

Deletes all the derived quantities.

Reimplemented in Lorene::AltBH_QI, Lorene::Boson_star, Lorene::Compobj_QI, Lorene::Kerr_QI, Lorene::ScalarBH, and Lorene::Star_QI.

Definition at line 155 of file compobj.C.

References p_adm_mass, and set_der_0x0().

◆ extrinsic_curvature()

void Lorene::Compobj::extrinsic_curvature ( )
virtual

Computation of the extrinsic curvature.

Reimplemented in Lorene::AltBH_QI, and Lorene::Compobj_QI.

Definition at line 290 of file compobj.C.

References beta, Lorene::Tensor::derive_cov(), gamma, kk, and nn.

◆ get_beta()

const Vector & Lorene::Compobj::get_beta ( ) const
inline

Returns the shift vector $\beta^i$.

Definition at line 213 of file compobj.h.

References beta.

◆ get_ener_euler()

const Scalar & Lorene::Compobj::get_ener_euler ( ) const
inline

Returns the total energy density E in the Eulerian frame.

Definition at line 219 of file compobj.h.

References ener_euler.

◆ get_gamma()

const Metric & Lorene::Compobj::get_gamma ( ) const
inline

Returns the 3-metric $\gamma_{ij}$.

Definition at line 216 of file compobj.h.

References gamma.

◆ get_kk()

const Sym_tensor & Lorene::Compobj::get_kk ( ) const
inline

Returns the extrinsic curvature tensor $K_{ij}$.

Definition at line 228 of file compobj.h.

References kk.

◆ get_mom_euler()

const Vector & Lorene::Compobj::get_mom_euler ( ) const
inline

Returns the total 3-momentum density $P^i$ in the Eulerian frame.

Definition at line 222 of file compobj.h.

References mom_euler.

◆ get_mp()

const Map & Lorene::Compobj::get_mp ( ) const
inline

Returns the mapping.

Definition at line 207 of file compobj.h.

References Lorene::Map(), and mp.

◆ get_nn()

const Scalar & Lorene::Compobj::get_nn ( ) const
inline

Returns the lapse function N .

Definition at line 210 of file compobj.h.

References nn.

◆ get_stress_euler()

const Sym_tensor & Lorene::Compobj::get_stress_euler ( ) const
inline

Returns the stress tensor $S_{ij}$ with respect to the Eulerian observer.

Definition at line 225 of file compobj.h.

References stress_euler.

◆ gyoto_data()

void Lorene::Compobj::gyoto_data ( const char * file_name) const

Save in a file for GYOTO.

Definition at line 210 of file compobj.C.

References beta, Lorene::fwrite_be(), gamma, kk, mp, and nn.

◆ operator=()

void Lorene::Compobj::operator= ( const Compobj & co)

Assignment to another Compobj.

Definition at line 175 of file compobj.C.

References beta, Compobj(), del_deriv(), ener_euler, gamma, kk, mom_euler, mp, nn, and stress_euler.

◆ operator>>()

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

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

Reimplemented in Lorene::AltBH_QI, Lorene::Boson_star, Lorene::Compobj_QI, Lorene::HiggsMonopole, Lorene::Kerr_QI, Lorene::ScalarBH, and Lorene::Star_QI.

Definition at line 239 of file compobj.C.

References ener_euler, gamma, kk, mp, nn, and stress_euler.

◆ sauve()

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

Save in a file.

Reimplemented in Lorene::AltBH_QI, Lorene::Boson_star, Lorene::Compobj_QI, Lorene::Kerr_QI, Lorene::ScalarBH, and Lorene::Star_QI.

Definition at line 196 of file compobj.C.

References beta, ener_euler, gamma, mom_euler, nn, and stress_euler.

◆ set_der_0x0()

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

Sets to 0x0 all the pointers on derived quantities.

Definition at line 163 of file compobj.C.

References p_adm_mass.

◆ set_mp()

Map & Lorene::Compobj::set_mp ( )
inline

Read/write of the mapping.

Definition at line 200 of file compobj.h.

References Lorene::Map(), and mp.

◆ operator<<

ostream & operator<< ( ostream & ost,
const Compobj & co )
friend

Display.

Definition at line 233 of file compobj.C.

References Compobj(), and operator<<.

Member Data Documentation

◆ beta

Vector Lorene::Compobj::beta
protected

Shift vector $\beta^i$.

Definition at line 138 of file compobj.h.

◆ ener_euler

Scalar Lorene::Compobj::ener_euler
protected

Total energy density E in the Eulerian frame.

Definition at line 144 of file compobj.h.

◆ gamma

Metric Lorene::Compobj::gamma
protected

3-metric $\gamma_{ij}$

Definition at line 141 of file compobj.h.

◆ kk

Sym_tensor Lorene::Compobj::kk
protected

Extrinsic curvature tensor $K_{ij}$.

Definition at line 153 of file compobj.h.

◆ mom_euler

Vector Lorene::Compobj::mom_euler
protected

Total 3-momentum density $P^i$ in the Eulerian frame.

Definition at line 147 of file compobj.h.

◆ mp

Map& Lorene::Compobj::mp
protected

Mapping describing the coordinate system (r,theta,phi).

Definition at line 132 of file compobj.h.

◆ nn

Scalar Lorene::Compobj::nn
protected

Lapse function N .

Definition at line 135 of file compobj.h.

◆ p_adm_mass

double* Lorene::Compobj::p_adm_mass
mutableprotected

ADM mass.

Definition at line 158 of file compobj.h.

◆ stress_euler

Sym_tensor Lorene::Compobj::stress_euler
protected

Stress tensor $S_{ij}$ with respect to the Eulerian observer.

Definition at line 150 of file compobj.h.


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