26 #ifndef __EXCISIONHOR_H_ 27 #define __EXCISIONHOR_H_ 189 const Scalar&
get_BC_bmN(
int choice_bmN,
double value = 1.)
const ;
192 const Scalar&
get_BC_bpN(
int choice_bpN,
double c_bpn_lap = 1.,
double c_bpn_fin = 1., Scalar *bpN_fin = 0x0)
const ;
199 virtual void sauve(FILE *)
const ;
202 friend ostream&
operator<<(ostream& ,
const Spheroid& ) ;
206 ostream& operator<<(ostream& ,
const Spheroid& ) ;
virtual void sauve(FILE *) const
Save in a file.
const Scalar & get_BC_conf_fact() const
Source of Neumann BC on , derived from the vanishing expansion.
double no_of_steps
The internal number of timesteps for one iteration.
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.
Metric for tensor calculation.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
const Spheroid & get_sph() const
Returns the spheroid.
const Sym_tensor & get_Tij() const
Returns the value of the impulsion-energy tensor.
Scalar & set_conf_fact()
Sets the value of the conformal factor.
virtual void del_deriv() const
Deletes all the derived quantities.
Tensor field of valence 0 (or component of a tensorial field).
Vector * p_get_BC_shift
Source of Dirichlet BC for the shift, issued from BC on bpN and a gauge condition on the tangential s...
Scalar & set_lapse()
Sets the lapse function.
const Scalar & get_lapse() const
Returns the lapse function.
const Sym_tensor & get_Kij() const
returns the 3-d extrinsic curvature
const Metric & get_gamij() const
Returns the symmetric tensor .
Vector shift
The Shift 3-vector on the slice.
Tensor field of valence 1.
const Scalar & get_conf_fact() const
Returns the conformal factor associated with the surface.
Scalar * p_get_BC_bmN
Source of Dirichlet BC for (b-N).
Surface where boundary conditions for quantities in the bulk will be calculated It relies on geometri...
friend ostream & operator<<(ostream &, const Spheroid &)
Display.
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Vector & set_shift()
Sets the shift vector field.
Sym_tensor & set_Kij()
Sets the extrinsic curvature.
Sym_tensor & set_Tij()
Sets the value of the impulsion-energy tensor.
const Vector & get_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 s...
Metric gamij
The 3-d metric on the slice.
void operator=(const Excision_hor &)
Assignment to another Excision_hor.
const Scalar & get_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 va...
Scalar * p_get_BC_conf_fact
Source of Neumann BC on , derived from the vanishing expansion.
Spheroid sph
The associated Spheroid object.
virtual ~Excision_hor()
Destructor.
Scalar * p_get_BC_bpN
Arbitrary source of Dirichlet BC for (b+N), case 0: based on a parabolic driver towards a constant va...
Metric & set_gamij()
Sets the 3d metric of the TimeSlice.
Sym_tensor Kij
The 3-d extrinsic curvature on the slice.
Sym_tensor Tij
Value of the impulsion-energy tensor on the spheroid.
const Vector & get_shift() const
Returns the shift vector field.
Scalar lapse
The lapse defined on the 3 slice.
double delta_t
The time step for evolution in parabolic drivers.
double get_delta_t() const
Returns the timestep used for evolution.
const Scalar & get_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 ...
Scalar conf_fact
The value of the conformal factor on the 3-slice.
Class intended to describe valence-2 symmetric tensors.
double get_no_of_steps() const
Returns the internal number of timesteps for one iteration.