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. | |
| Excision_hor (const Excision_hor &) | |
| Copy constructor. | |
| Excision_hor (FILE *) | |
Constructor from a file (see sauve(FILE*) ). | |
| virtual | ~Excision_hor () |
| Destructor. | |
| void | operator= (const Excision_hor &) |
| Assignment to another Excision_hor. | |
| const Spheroid & | get_sph () const |
| Returns the spheroid. | |
| const Scalar & | get_conf_fact () const |
| Returns the conformal factor associated with the surface. | |
| const Scalar & | get_lapse () const |
| Returns the lapse function. | |
| const Vector & | get_shift () const |
| Returns the shift vector field. | |
| const Metric & | get_gamij () const |
Returns the symmetric tensor . | |
| const Sym_tensor & | get_Kij () const |
returns the 3-d extrinsic curvature | |
| double | get_delta_t () const |
| Returns the timestep used for evolution. | |
| double | get_no_of_steps () const |
| Returns the internal number of timesteps for one iteration. | |
| 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. | |
| Scalar & | set_lapse () |
| Sets the lapse function. | |
| Vector & | set_shift () |
| Sets the shift vector field. | |
| Metric & | set_gamij () |
| Sets the 3d metric of the TimeSlice. | |
| Sym_tensor & | set_Kij () |
| Sets the extrinsic curvature. | |
| Sym_tensor & | set_Tij () |
| Sets the value of the impulsion-energy tensor. | |
| double & | set_delta_t () |
| double & | set_no_of_steps () |
| const Scalar & | get_BC_conf_fact () const |
Source of Neumann BC on , derived from the vanishing expansion. | |
| 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 of projected Einstein Equations. | |
| 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 value. | |
| 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 shift (based on a parabolic driver). | |
| virtual void | sauve (FILE *) const |
| Save in a file. | |
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. | |
Protected Attributes | |
| Spheroid | sph |
| The associated Spheroid object. | |
| Scalar | conf_fact |
| The value of the conformal factor on the 3-slice. | |
| Scalar | lapse |
| The lapse defined on the 3 slice. | |
| Vector | shift |
| The Shift 3-vector on the slice. | |
| Metric | gamij |
| The 3-d metric on the slice. | |
| Sym_tensor | Kij |
| The 3-d extrinsic curvature on the slice. | |
| double | delta_t |
| The time step for evolution in parabolic drivers. | |
| double | no_of_steps |
| The internal number of timesteps for one iteration. | |
| Sym_tensor | Tij |
| Value of the impulsion-energy tensor on the spheroid. | |
| Scalar * | p_get_BC_conf_fact |
Source of Neumann BC on , derived from the vanishing expansion. | |
| Scalar * | p_get_BC_bmN |
| Source of Dirichlet BC for (b-N). | |
| Scalar * | p_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. . | |
| Vector * | p_get_BC_shift |
| Source of Dirichlet BC for the shift, issued from BC on bpN and a gauge condition on the tangential shift. | |
Friends | |
| ostream & | operator<< (ostream &, const Spheroid &) |
| Display. | |
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 42 of file excision_hor.h.
| 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.
| 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 45 of file excision_hor.C.
References set_der_0x0().
| Excision_hor::Excision_hor | ( | const Excision_hor & | exc_in | ) |
| Excision_hor::Excision_hor | ( | FILE * | ) |
Constructor from a file (see sauve(FILE*) ).
| Excision_hor::~Excision_hor | ( | ) | [virtual] |
| void Excision_hor::del_deriv | ( | ) | const [protected, virtual] |
Deletes all the derived quantities.
Definition at line 93 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().
| const Scalar & Excision_hor::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 of projected Einstein Equations.
Definition at line 146 of file excision_hor.C.
References Scalar::allocate_all(), Metric_flat::con(), Metric::con(), contract(), Spheroid::delta(), Spheroid::derive_cov2d(), Spheroid::derive_cov2dflat(), Map::flat_met_spher(), gamij, Spheroid::get_hsurf(), Spheroid::get_ll(), Map::get_mg(), Tensor::get_mp(), Spheroid::get_qab(), Spheroid::get_ricci(), lapse, max(), maxabs(), p_get_BC_bmN, Scalar::poisson_angu(), Metric::radial_vect(), Tensor::set(), Scalar::set_spectral_va(), Spheroid::shear(), shift, sph, Spheroid::sqrt_q(), Tensor::std_spectral_base(), Spheroid::theta_minus(), Tij, Tensor::trace(), Tensor::up_down(), Scalar::val_grid_point(), and Valeur::ylm().
| const Scalar & Excision_hor::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 value.
Case 1: Source of Dirichlet BC for (b+N), from a component of projected Einstein Equations.
Definition at line 297 of file excision_hor.C.
References contract(), delta_t, Spheroid::derive_cov2d(), gamij, get_BC_bmN(), Spheroid::get_hsurf(), Spheroid::get_ll(), Map::get_mg(), Tensor::get_mp(), Spheroid::get_qab(), Spheroid::get_ricci(), Scalar::lapang(), lapse, p_get_BC_bpN, Metric::radial_vect(), Scalar::set_spectral_va(), Spheroid::shear(), shift, sph, Scalar::std_spectral_base(), Tij, Tensor::trace(), Tensor::up(), Tensor::up_down(), Scalar::val_grid_point(), and Valeur::ylm().
| const Scalar & Excision_hor::get_BC_conf_fact | ( | ) | const |
Source of Neumann BC on
, derived from the vanishing expansion.
Definition at line 119 of file excision_hor.C.
References conf_fact, contract(), Metric::cov(), Scalar::derive_cov(), Vector::divergence(), Scalar::dsdr(), gamij, Kij, p_get_BC_conf_fact, pow(), Metric::radial_vect(), Scalar::set_spectral_va(), Scalar::std_spectral_base(), Tensor::std_spectral_base(), and Valeur::ylm().
| const Vector & Excision_hor::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 shift (based on a parabolic driver).
Definition at line 420 of file excision_hor.C.
References Scalar::allocate_all(), Tensor::annule_domain(), Metric::con(), contract(), Metric::cov(), delta_t, Tensor::derive_con(), Tensor_sym::derive_cov(), Tensor::derive_cov(), Tensor::down(), gamij, Map::get_mg(), Tensor::get_mp(), Spheroid::get_ricci(), lapse, p_get_BC_bmN, p_get_BC_bpN, p_get_BC_shift, Metric::radial_vect(), shift, sph, Vector::std_spectral_base(), Tensor::std_spectral_base(), Tensor::up_down(), and Scalar::val_grid_point().
| const Scalar& Excision_hor::get_conf_fact | ( | ) | const [inline] |
Returns the conformal factor associated with the surface.
Definition at line 134 of file excision_hor.h.
References conf_fact.
| double Excision_hor::get_delta_t | ( | ) | const [inline] |
Returns the timestep used for evolution.
Definition at line 149 of file excision_hor.h.
References delta_t.
| const Metric& Excision_hor::get_gamij | ( | ) | const [inline] |
| const Sym_tensor& Excision_hor::get_Kij | ( | ) | const [inline] |
| const Scalar& Excision_hor::get_lapse | ( | ) | const [inline] |
| double Excision_hor::get_no_of_steps | ( | ) | const [inline] |
Returns the internal number of timesteps for one iteration.
Definition at line 152 of file excision_hor.h.
References no_of_steps.
| const Vector& Excision_hor::get_shift | ( | ) | const [inline] |
| const Spheroid& Excision_hor::get_sph | ( | ) | const [inline] |
| const Sym_tensor& Excision_hor::get_Tij | ( | ) | const [inline] |
Returns the value of the impulsion-energy tensor.
Definition at line 155 of file excision_hor.h.
References Tij.
| void Excision_hor::operator= | ( | const Excision_hor & | ) |
Assignment to another Excision_hor.
| void Excision_hor::sauve | ( | FILE * | ) | const [virtual] |
Save in a file.
Definition at line 531 of file excision_hor.C.
| Scalar& Excision_hor::set_conf_fact | ( | ) | [inline] |
Sets the value of the conformal factor.
Definition at line 158 of file excision_hor.h.
References conf_fact, and del_deriv().
| void Excision_hor::set_der_0x0 | ( | ) | const [protected] |
Sets to 0x0 all the pointers on derived quantities.
Definition at line 102 of file excision_hor.C.
References p_get_BC_bmN, p_get_BC_bpN, p_get_BC_conf_fact, and p_get_BC_shift.
| Metric& Excision_hor::set_gamij | ( | ) | [inline] |
Sets the 3d metric of the TimeSlice.
Definition at line 167 of file excision_hor.h.
References del_deriv(), and gamij.
| Sym_tensor& Excision_hor::set_Kij | ( | ) | [inline] |
Sets the extrinsic curvature.
Definition at line 170 of file excision_hor.h.
References del_deriv(), and Kij.
| Scalar& Excision_hor::set_lapse | ( | ) | [inline] |
Sets the lapse function.
Definition at line 161 of file excision_hor.h.
References del_deriv(), and lapse.
| Vector& Excision_hor::set_shift | ( | ) | [inline] |
Sets the shift vector field.
Definition at line 164 of file excision_hor.h.
References del_deriv(), and shift.
| Sym_tensor& Excision_hor::set_Tij | ( | ) | [inline] |
Sets the value of the impulsion-energy tensor.
Definition at line 173 of file excision_hor.h.
References del_deriv(), and Tij.
| ostream& operator<< | ( | ostream & | , | |
| const Spheroid & | ||||
| ) | [friend] |
Display.
Scalar Excision_hor::conf_fact [protected] |
The value of the conformal factor on the 3-slice.
Definition at line 52 of file excision_hor.h.
double Excision_hor::delta_t [protected] |
The time step for evolution in parabolic drivers.
Definition at line 67 of file excision_hor.h.
Metric Excision_hor::gamij [protected] |
The 3-d metric on the slice.
Definition at line 61 of file excision_hor.h.
Sym_tensor Excision_hor::Kij [protected] |
The 3-d extrinsic curvature on the slice.
Definition at line 64 of file excision_hor.h.
Scalar Excision_hor::lapse [protected] |
The lapse defined on the 3 slice.
Definition at line 55 of file excision_hor.h.
double Excision_hor::no_of_steps [protected] |
The internal number of timesteps for one iteration.
Definition at line 70 of file excision_hor.h.
Scalar* Excision_hor::p_get_BC_bmN [mutable, protected] |
Source of Dirichlet BC for (b-N).
Definition at line 80 of file excision_hor.h.
Scalar* Excision_hor::p_get_BC_bpN [mutable, protected] |
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 81 of file excision_hor.h.
Scalar* Excision_hor::p_get_BC_conf_fact [mutable, protected] |
Source of Neumann BC on
, derived from the vanishing expansion.
Definition at line 79 of file excision_hor.h.
Vector* Excision_hor::p_get_BC_shift [mutable, protected] |
Source of Dirichlet BC for the shift, issued from BC on bpN and a gauge condition on the tangential shift.
Definition at line 82 of file excision_hor.h.
Vector Excision_hor::shift [protected] |
The Shift 3-vector on the slice.
Definition at line 58 of file excision_hor.h.
Spheroid Excision_hor::sph [protected] |
The associated Spheroid object.
Definition at line 49 of file excision_hor.h.
Sym_tensor Excision_hor::Tij [protected] |
Value of the impulsion-energy tensor on the spheroid.
Definition at line 73 of file excision_hor.h.
1.6.1