LORENE
Lorene::Excision_hor Class Reference

Surface where boundary conditions for quantities in the bulk will be calculated It relies on geometrical properties of the associated Spheroid() (*** WARNING! under development***) More...

#include <excision_hor.h>

Public Member Functions

 Excision_hor (const Scalar &h_in, const Metric &gij, const Sym_tensor &Kij2, const Scalar &ppsi, const Scalar &nn, const Vector &beta, const Sym_tensor &Tij2, double timestep, int int_nos=1)
 Constructor of an excision surface embedded in a 3-slice (Time_slice ) of 3+1 formalism. More...
 
 Excision_hor (const Excision_hor &)
 Copy constructor. More...
 
 Excision_hor (FILE *)
 Constructor from a file (see sauve(FILE*) ) More...
 
virtual ~Excision_hor ()
 Destructor. More...
 
void operator= (const Excision_hor &)
 Assignment to another Excision_hor. More...
 
const Spheroidget_sph () const
 Returns the spheroid. More...
 
const Scalarget_conf_fact () const
 Returns the conformal factor associated with the surface. More...
 
const Scalarget_lapse () const
 Returns the lapse function. More...
 
const Vectorget_shift () const
 Returns the shift vector field. More...
 
const Metricget_gamij () const
 Returns the symmetric tensor $ gamij $. More...
 
const Sym_tensorget_Kij () const
 returns the 3-d extrinsic curvature $ K_{ij}$ More...
 
double get_delta_t () const
 Returns the timestep used for evolution. More...
 
double get_no_of_steps () const
 Returns the internal number of timesteps for one iteration. More...
 
const Sym_tensorget_Tij () const
 Returns the value of the impulsion-energy tensor. More...
 
Scalarset_conf_fact ()
 Sets the value of the conformal factor. More...
 
Scalarset_lapse ()
 Sets the lapse function. More...
 
Vectorset_shift ()
 Sets the shift vector field. More...
 
Metricset_gamij ()
 Sets the 3d metric of the TimeSlice. More...
 
Sym_tensorset_Kij ()
 Sets the extrinsic curvature. More...
 
Sym_tensorset_Tij ()
 Sets the value of the impulsion-energy tensor. More...
 
double & set_delta_t ()
 
double & set_no_of_steps ()
 
const Scalarget_BC_conf_fact () const
 Source of Neumann BC on $ \psi $, derived from the vanishing expansion. More...
 
const Scalarget_BC_bmN (int choice_bmN, double value=1.) const
 Source of Dirichlet BC for (b-N): case 0: based on an entropy prescription, case 1: from a component of projected Einstein Equations. More...
 
const Scalarget_BC_bpN (int choice_bpN, double c_bpn_lap=1., double c_bpn_fin=1., Scalar *bpN_fin=0x0) const
 Case 0: Arbitrary source of Dirichlet BC for (b+N), based on a parabolic driver towards a constant value. More...
 
const Vectorget_BC_shift (double c_V_lap) const
 Source of Dirichlet BC for the shift, issued from BC on bpN and a gauge condition on the tangential shift (based on a parabolic driver). More...
 
virtual void sauve (FILE *) const
 Save in a file. 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...
 

Protected Attributes

Spheroid sph
 The associated Spheroid object. More...
 
Scalar conf_fact
 The value of the conformal factor on the 3-slice. More...
 
Scalar lapse
 The lapse defined on the 3 slice. More...
 
Vector shift
 The Shift 3-vector on the slice. More...
 
Metric gamij
 The 3-d metric on the slice. More...
 
Sym_tensor Kij
 The 3-d extrinsic curvature on the slice. More...
 
double delta_t
 The time step for evolution in parabolic drivers. More...
 
double no_of_steps
 The internal number of timesteps for one iteration. More...
 
Sym_tensor Tij
 Value of the impulsion-energy tensor on the spheroid. More...
 
Scalarp_get_BC_conf_fact
 Source of Neumann BC on $ \psi $, derived from the vanishing expansion. More...
 
Scalarp_get_BC_bmN
 Source of Dirichlet BC for (b-N). More...
 
Scalarp_get_BC_bpN
 Arbitrary source of Dirichlet BC for (b+N), case 0: based on a parabolic driver towards a constant value, case 1:from a component of projected Einstein Equations. . More...
 
Vectorp_get_BC_shift
 Source of Dirichlet BC for the shift, issued from BC on bpN and a gauge condition on the tangential shift. More...
 

Friends

ostream & operator<< (ostream &, const Spheroid &)
 Display. More...
 

Detailed Description

Surface where boundary conditions for quantities in the bulk will be calculated It relies on geometrical properties of the associated Spheroid() (*** WARNING! under development***)

Definition at line 43 of file excision_hor.h.

Constructor & Destructor Documentation

◆ Excision_hor() [1/3]

Lorene::Excision_hor::Excision_hor ( const Scalar h_in,
const Metric gij,
const Sym_tensor Kij2,
const Scalar ppsi,
const Scalar nn,
const Vector beta,
const Sym_tensor Tij2,
double  timestep,
int  int_nos = 1 
)

Constructor of an excision surface embedded in a 3-slice (Time_slice ) of 3+1 formalism.

This is done from the Time_slice data.

Parameters
h_in: the location of the surface r = h_in (WARNING:must be defined on a mono-domain angular grid)
gij: the 3-metric on the 3-slice
Kij: the extrinsic curvature of the 3-slice (covariant representation)
timestep: time interval associated with the parabolic-driven boundary conditions.
int_nos: Number of iterations to be done during timestep.
Tij: Value of the impulsion-energy tensor on the spheroid

Definition at line 46 of file excision_hor.C.

References set_der_0x0().

◆ Excision_hor() [2/3]

Lorene::Excision_hor::Excision_hor ( const Excision_hor exc_in)

Copy constructor.

Definition at line 68 of file excision_hor.C.

References set_der_0x0().

◆ Excision_hor() [3/3]

Lorene::Excision_hor::Excision_hor ( FILE *  )

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

◆ ~Excision_hor()

Lorene::Excision_hor::~Excision_hor ( )
virtual

Destructor.

Definition at line 86 of file excision_hor.C.

References del_deriv().

Member Function Documentation

◆ del_deriv()

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

Deletes all the derived quantities.

Definition at line 94 of file excision_hor.C.

References p_get_BC_bmN, p_get_BC_bpN, p_get_BC_conf_fact, p_get_BC_shift, and set_der_0x0().

◆ get_BC_bmN()

◆ get_BC_bpN()

const Scalar & Lorene::Excision_hor::get_BC_bpN ( int  choice_bpN,
double  c_bpn_lap = 1.,
double  c_bpn_fin = 1.,
Scalar bpN_fin = 0x0 
) const

◆ get_BC_conf_fact()

◆ get_BC_shift()

◆ get_conf_fact()

const Scalar& Lorene::Excision_hor::get_conf_fact ( ) const
inline

Returns the conformal factor associated with the surface.

Definition at line 135 of file excision_hor.h.

References conf_fact.

◆ get_delta_t()

double Lorene::Excision_hor::get_delta_t ( ) const
inline

Returns the timestep used for evolution.

Definition at line 150 of file excision_hor.h.

References delta_t.

◆ get_gamij()

const Metric& Lorene::Excision_hor::get_gamij ( ) const
inline

Returns the symmetric tensor $ gamij $.

Definition at line 144 of file excision_hor.h.

References gamij.

◆ get_Kij()

const Sym_tensor& Lorene::Excision_hor::get_Kij ( ) const
inline

returns the 3-d extrinsic curvature $ K_{ij}$

Definition at line 147 of file excision_hor.h.

References Kij.

◆ get_lapse()

const Scalar& Lorene::Excision_hor::get_lapse ( ) const
inline

Returns the lapse function.

Definition at line 138 of file excision_hor.h.

References lapse.

◆ get_no_of_steps()

double Lorene::Excision_hor::get_no_of_steps ( ) const
inline

Returns the internal number of timesteps for one iteration.

Definition at line 153 of file excision_hor.h.

References no_of_steps.

◆ get_shift()

const Vector& Lorene::Excision_hor::get_shift ( ) const
inline

Returns the shift vector field.

Definition at line 141 of file excision_hor.h.

References shift.

◆ get_sph()

const Spheroid& Lorene::Excision_hor::get_sph ( ) const
inline

Returns the spheroid.

Definition at line 132 of file excision_hor.h.

References sph.

◆ get_Tij()

const Sym_tensor& Lorene::Excision_hor::get_Tij ( ) const
inline

Returns the value of the impulsion-energy tensor.

Definition at line 156 of file excision_hor.h.

References Tij.

◆ operator=()

void Lorene::Excision_hor::operator= ( const Excision_hor )

Assignment to another Excision_hor.

◆ sauve()

void Lorene::Excision_hor::sauve ( FILE *  ) const
virtual

Save in a file.

Definition at line 540 of file excision_hor.C.

◆ set_conf_fact()

Scalar& Lorene::Excision_hor::set_conf_fact ( )
inline

Sets the value of the conformal factor.

Definition at line 159 of file excision_hor.h.

References conf_fact, and del_deriv().

◆ set_der_0x0()

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

Sets to 0x0 all the pointers on derived quantities.

Definition at line 103 of file excision_hor.C.

References p_get_BC_bmN, p_get_BC_bpN, p_get_BC_conf_fact, and p_get_BC_shift.

◆ set_gamij()

Metric& Lorene::Excision_hor::set_gamij ( )
inline

Sets the 3d metric of the TimeSlice.

Definition at line 168 of file excision_hor.h.

References del_deriv(), and gamij.

◆ set_Kij()

Sym_tensor& Lorene::Excision_hor::set_Kij ( )
inline

Sets the extrinsic curvature.

Definition at line 171 of file excision_hor.h.

References del_deriv(), and Kij.

◆ set_lapse()

Scalar& Lorene::Excision_hor::set_lapse ( )
inline

Sets the lapse function.

Definition at line 162 of file excision_hor.h.

References del_deriv(), and lapse.

◆ set_shift()

Vector& Lorene::Excision_hor::set_shift ( )
inline

Sets the shift vector field.

Definition at line 165 of file excision_hor.h.

References del_deriv(), and shift.

◆ set_Tij()

Sym_tensor& Lorene::Excision_hor::set_Tij ( )
inline

Sets the value of the impulsion-energy tensor.

Definition at line 174 of file excision_hor.h.

References del_deriv(), and Tij.

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  ,
const Spheroid  
)
friend

Display.

Member Data Documentation

◆ conf_fact

Scalar Lorene::Excision_hor::conf_fact
protected

The value of the conformal factor on the 3-slice.

Definition at line 53 of file excision_hor.h.

◆ delta_t

double Lorene::Excision_hor::delta_t
protected

The time step for evolution in parabolic drivers.

Definition at line 68 of file excision_hor.h.

◆ gamij

Metric Lorene::Excision_hor::gamij
protected

The 3-d metric on the slice.

Definition at line 62 of file excision_hor.h.

◆ Kij

Sym_tensor Lorene::Excision_hor::Kij
protected

The 3-d extrinsic curvature on the slice.

Definition at line 65 of file excision_hor.h.

◆ lapse

Scalar Lorene::Excision_hor::lapse
protected

The lapse defined on the 3 slice.

Definition at line 56 of file excision_hor.h.

◆ no_of_steps

double Lorene::Excision_hor::no_of_steps
protected

The internal number of timesteps for one iteration.

Definition at line 71 of file excision_hor.h.

◆ p_get_BC_bmN

Scalar* Lorene::Excision_hor::p_get_BC_bmN
mutableprotected

Source of Dirichlet BC for (b-N).

Definition at line 81 of file excision_hor.h.

◆ p_get_BC_bpN

Scalar* Lorene::Excision_hor::p_get_BC_bpN
mutableprotected

Arbitrary source of Dirichlet BC for (b+N), case 0: based on a parabolic driver towards a constant value, case 1:from a component of projected Einstein Equations. .

Definition at line 82 of file excision_hor.h.

◆ p_get_BC_conf_fact

Scalar* Lorene::Excision_hor::p_get_BC_conf_fact
mutableprotected

Source of Neumann BC on $ \psi $, derived from the vanishing expansion.

Definition at line 80 of file excision_hor.h.

◆ p_get_BC_shift

Vector* Lorene::Excision_hor::p_get_BC_shift
mutableprotected

Source of Dirichlet BC for the shift, issued from BC on bpN and a gauge condition on the tangential shift.

Definition at line 83 of file excision_hor.h.

◆ shift

Vector Lorene::Excision_hor::shift
protected

The Shift 3-vector on the slice.

Definition at line 59 of file excision_hor.h.

◆ sph

Spheroid Lorene::Excision_hor::sph
protected

The associated Spheroid object.

Definition at line 50 of file excision_hor.h.

◆ Tij

Sym_tensor Lorene::Excision_hor::Tij
protected

Value of the impulsion-energy tensor on the spheroid.

Definition at line 74 of file excision_hor.h.


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