26 #ifndef __EXCISIONSURF_H_ 27 #define __EXCISIONSURF_H_ 232 const Scalar& get_BC_lapse_1(
double value)
const ;
234 const Vector& get_BC_shift_1(
double Omega)
const ;
236 const Scalar& get_BC_Npsi_1(
double value)
const ;
252 const Vector&
get_BC_shift_2(
double c_bb_lap,
double c_bb_fin,
double c_V_lap,
double epsilon)
const ;
268 virtual void sauve(FILE *)
const ;
275 ostream& operator<<(ostream& ,
const Spheroid& ) ;
const Scalar & get_BC_Npsi_5(double Kappa) const
Source for a Neumann BC on (N*Psi1), fixing a constant surface gravity in space and time...
Scalar * p_get_BC_lapse_4
Source of Dirichlet condtion on , based on einstein equations (conservation of isotropic gauge) ...
Metric & set_gamij()
Sets the 3d metric of the TimeSlice.
Metric for tensor calculation.
void set_expa_hyperb(double alph0, double beta0, double gamma0)
Sets a new value for expansion rescaled over lapse (and its derivative), obtained by hyperbolic evolu...
const Scalar & get_BC_conf_fact_2(double c_psi_lap, double c_psi_fin, Scalar &expa_fin) const
Source for the Dirichlet BC on the conformal factor, based on a parabolic driver for the conformal fa...
double delta_t
The time step for evolution in parabolic drivers.
Sym_tensor Kij
The 3-d extrinsic curvature on the slice.
Scalar * p_get_BC_Npsi_2
Source of Dirichlet boundary condition on .
Scalar * p_get_BC_Npsi_4
Source of Dirichlet boundary condition on .
const Scalar & derive_t_expa(Scalar &Ee, Vector &Jj, Sym_tensor &Sij) const
Forms the prospective time derivative for the expansion using projected Einstein equations. Does NOT modify the member dt_expa: do it by hand!
const Scalar & get_expa() const
Returns the assumed expansion associated to the excised surface at t.
Scalar expa
The 2d expansion, directly evolved from the initial excision with Einstein Equations.
const Spheroid & get_sph() const
Returns the spheroid.
Vector * p_get_BC_shift_4
Source of Dirichlet BC for the shift vector , partly from projection of Einstein Equations.
const Vector & get_BC_shift_3(Scalar &dtpsi, double c_V_lap, double epsilon) const
Source for a Dirichlet BC on the shift, based on a Parabolic driver; Radial part is dealt with using ...
Scalar * p_get_BC_conf_fact_2
Source of Neumann boundary condition on ,.
const Scalar & get_BC_conf_fact_3(double c_theta_lap, double c_theta_fin, Scalar &expa_fin) const
Source for the Neumann BC on the conformal factor, based on a parabolic driver for the expansio...
const Scalar & get_BC_Npsi_4(double Kappa) const
Source for a Dirichlet BC on (N*Psi1), fixing a constant surface gravity in space and time...
Scalar & set_expa()
Sets the expansion function on the surface at time t (considering to protect this function) ...
Tensor field of valence 0 (or component of a tensorial field).
void operator=(const Excision_surf &)
Assignment to another Excision_surf.
Spheroid sph
The associated Spheroid object.
Tensor field of valence 1.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Scalar & set_dt_expa()
Sets the time derivative of the expansion function on the surface at time t (considering to protect t...
Sym_tensor & set_Kij()
Sets the extrinsic curvature.
Scalar * p_get_BC_lapse_3
Source of Dirichlet condtion on , based on einstein equations.
Scalar * p_get_BC_Npsi_5
Source of Neumann boundary condition on .
Scalar lapse
The lapse defined on the 3 slice.
const Sym_tensor & get_Kij() const
returns the 3-d extrinsic curvature
Excision_surf(const Scalar &h_in, const Metric &gij, const Sym_tensor &Kij2, const Scalar &ppsi, const Scalar &nn, const Vector &beta, double timestep, int int_nos)
Constructor of an excision surface embedded in a 3-slice (Time_slice ) of 3+1 formalism.
const Metric & get_gamij() const
Returns the symmetric tensor .
const Scalar & get_dt_expa() const
Returns the assumed time derivative of the expansion at t, evolved by functions of this class;...
void set_expa_parab(double c_theta_lap, double c_theta_fin, Scalar &expa_fin)
Sets a new value for expansion rescaled over lapse (and its derivative), obtained by parabolic evolut...
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Scalar & set_lapse()
Sets the lapse function.
Scalar * p_get_BC_conf_fact_4
Source of Birichlet boundary condition on ,.
Scalar * p_get_BC_conf_fact_3
Source of Neumann boundary condition on ,.
virtual void del_deriv() const
Deletes all the derived quantities.
double get_no_of_steps() const
Returns the internal number of timesteps for one iteration.
void get_evol_params_from_ID(double alpha, double beta, double gamma, Scalar &Ee, Vector &Jj, Sym_tensor &Ss)
Computes the parameters for the hyperbolic evolution in set_expa_hyperb(), so that the expansion has ...
Scalar * p_get_BC_Npsi_3
Source of Dirichlet boundary condition on .
Vector & set_shift()
Sets the shift vector field.
const Vector & get_shift() const
Returns the shift vector field.
Scalar * p_derive_t_expa
Computation of an updated expansion scalar.
const Scalar & get_conf_fact() const
Returns the conformal factor associated with the surface.
virtual void sauve(FILE *) const
Save in a file.
const Scalar & get_BC_Npsi_2(double value, double c_npsi_lap, double c_npsi_fin) const
Source for the Dirichlet BC on (N*Psi1), based on a parabolic driver.
const Scalar & get_lapse() const
Returns the lapse function.
Vector * p_get_BC_shift_3
Source of Dirichlet BC for the shift vector , partly derived from kinematical relation.
Spheroid & set_sph()
Sets a new spheroid from data.
double get_delta_t() const
Returns the timestep used for evolution.
const Scalar & get_BC_lapse_3(Scalar &dttheta, Scalar &Ee, Vector &Jj, Sym_tensor &Sij, bool sph_sym=true) const
Source for Dirichlet BC on the lapse, based on einstein equations.
Scalar * p_get_BC_lapse_1
Source of Dirichlet boundary condition of .
double no_of_steps
The internal number of timesteps for one iteration.
Surface where boundary conditions for quantities in the bulk will be calculated It relies on geometri...
virtual ~Excision_surf()
Destructor.
Scalar dt_expa
The time derivative of the expansion, derived from Einstein equations and arbitrary evolution...
Scalar conf_fact
The value of the conformal factor on the 3-slice.
const Scalar & get_BC_lapse_4(Scalar &old_nn, Vector &beta_point, Sym_tensor &strain_tens) const
Source for Dirichlet BC on the lapse, based on einstein equations (conservation of isotropic gauge) ...
const Scalar & get_BC_Npsi_3(double n_0, double beta) const
Source for the Dirichlet BC on (N*Psi1), with Kerr_Schild-like form for the lapse boundary...
const Vector & get_BC_shift_2(double c_bb_lap, double c_bb_fin, double c_V_lap, double epsilon) const
Source for a Dirichlet BC on the shift, based on a Parabolic driver; no assumptions are made except a...
Scalar * p_get_BC_lapse_2
Source of Dirichlet boundary condition of .
const Scalar & get_BC_lapse_2(double lapse_fin, double c_lapse_lap, double c_lapse_fi) const
Source for Dirichlet BC on the lapse, based on a parabolic driver towards arbitrary constant value...
const Scalar & get_BC_conf_fact_4() const
Source for the Dirchlet BC on the conformal factor, based on the consistency condition derived from t...
Scalar & set_conf_fact()
Sets the value of the conformal factor.
Vector * p_get_BC_shift_1
Source of Dirichlet BC for the shift vector .
Metric gamij
The 3-d metric on the slice.
Vector * p_get_BC_shift_2
Source of Dirichlet BC for the shift vector .
Scalar * p_get_BC_Npsi_1
Source of Neumann boundary condition on .
friend ostream & operator<<(ostream &, const Spheroid &)
Display.
const Vector & get_BC_shift_4(Scalar &dttheta, Scalar &Ee, Vector &Jj, Sym_tensor &Sij, double c_V_lap, double epsilon, bool sph_sym=true) const
Source for a Dirichlet BC on the shift, based on a Parabolic driver; Radial part is dealt with using ...
Vector shift
The Shift 3-vector on the slice.
Class intended to describe valence-2 symmetric tensors.
const Scalar & get_BC_conf_fact_1(bool isMOTS=false) const
Source for a Neumann BC on the conformal factor. If boolean isMOTS is false, it is based on expansion...
Scalar * p_get_BC_conf_fact_1
Source of Neumann boundary condition on ,.