LORENE
Lorene::Excision_surf 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_surf.h>

Public Member Functions

 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. More...
 
 Excision_surf (const Excision_surf &)
 Copy constructor. More...
 
 Excision_surf (FILE *)
 Constructor from a file (see sauve(FILE*) ) More...
 
virtual ~Excision_surf ()
 Destructor. More...
 
void operator= (const Excision_surf &)
 Assignment to another Excision_surf. More...
 
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 a C1 matching with initial data. More...
 
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 evolution. More...
 
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 evolution. 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 Scalarget_expa () const
 Returns the assumed expansion associated to the excised surface at t. More...
 
const Scalarget_dt_expa () const
 Returns the assumed time derivative of the expansion at t, evolved by functions of this class;. More...
 
Spheroidset_sph ()
 Sets a new spheroid from data. 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...
 
double set_delta_t ()
 
double set_no_of_steps ()
 
Scalarset_expa ()
 Sets the expansion function on the surface at time t (considering to protect this function) More...
 
Scalarset_dt_expa ()
 Sets the time derivative of the expansion function on the surface at time t (considering to protect this function) More...
 
const Scalarget_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 value of the spheroid or the value of exppa; it is based on zero expansion if isMOTS is true. More...
 
const Scalarget_BC_lapse_1 (double value) const
 
const Vectorget_BC_shift_1 (double Omega) const
 
const Scalarget_BC_Npsi_1 (double value) const
 
const Scalarget_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 factor. More...
 
const Scalarget_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. More...
 
const Scalarget_BC_conf_fact_4 () const
 Source for the Dirchlet BC on the conformal factor, based on the consistency condition derived from the trace. More...
 
const Scalarget_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. More...
 
const Scalarget_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. More...
 
const Scalarget_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) More...
 
const Scalarderive_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! More...
 
const Vectorget_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 global conformal Killing symmetry. More...
 
const Vectorget_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 a kinematical relation. More...
 
const Vectorget_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 projection of Einstein Equations. More...
 
const Scalarget_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. More...
 
const Scalarget_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. More...
 
const Scalarget_BC_Npsi_4 (double Kappa) const
 Source for a Dirichlet BC on (N*Psi1), fixing a constant surface gravity in space and time. More...
 
const Scalarget_BC_Npsi_5 (double Kappa) const
 Source for a Neumann BC on (N*Psi1), fixing a constant surface gravity in space and time. 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...
 
Scalar expa
 The 2d expansion, directly evolved from the initial excision with Einstein Equations. More...
 
Scalar dt_expa
 The time derivative of the expansion, derived from Einstein equations and arbitrary evolution. More...
 
Scalarp_get_BC_conf_fact_1
 Source of Neumann boundary condition on $ \psi $,. More...
 
Scalarp_get_BC_lapse_1
 Source of Dirichlet boundary condition of $ N $. More...
 
Vectorp_get_BC_shift_1
 Source of Dirichlet BC for the shift vector $ \beta^{i} $. More...
 
Scalarp_get_BC_Npsi_1
 Source of Neumann boundary condition on $ \psi $. More...
 
Scalarp_get_BC_conf_fact_2
 Source of Neumann boundary condition on $ \psi $,. More...
 
Scalarp_get_BC_conf_fact_3
 Source of Neumann boundary condition on $ \psi $,. More...
 
Scalarp_get_BC_conf_fact_4
 Source of Birichlet boundary condition on $ \psi $,. More...
 
Scalarp_get_BC_lapse_2
 Source of Dirichlet boundary condition of $ N $. More...
 
Scalarp_get_BC_lapse_3
 Source of Dirichlet condtion on $ N $, based on einstein equations. More...
 
Scalarp_get_BC_lapse_4
 Source of Dirichlet condtion on $ N $, based on einstein equations (conservation of isotropic gauge) More...
 
Scalarp_derive_t_expa
 Computation of an updated expansion scalar. More...
 
Vectorp_get_BC_shift_2
 Source of Dirichlet BC for the shift vector $ \beta^{i} $. More...
 
Vectorp_get_BC_shift_3
 Source of Dirichlet BC for the shift vector $ \beta^{i} $, partly derived from kinematical relation. More...
 
Vectorp_get_BC_shift_4
 Source of Dirichlet BC for the shift vector $ \beta^{i} $, partly from projection of Einstein Equations. More...
 
Scalarp_get_BC_Npsi_2
 Source of Dirichlet boundary condition on $ N \psi $. More...
 
Scalarp_get_BC_Npsi_3
 Source of Dirichlet boundary condition on $ N \psi $. More...
 
Scalarp_get_BC_Npsi_4
 Source of Dirichlet boundary condition on $ N \psi $. More...
 
Scalarp_get_BC_Npsi_5
 Source of Neumann boundary condition on $ N \psi $. 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_surf.h.

Constructor & Destructor Documentation

◆ Excision_surf() [1/3]

Lorene::Excision_surf::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 = 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 3slice
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.

Definition at line 64 of file excision_surf.C.

References dt_expa, set_der_0x0(), and Lorene::Scalar::set_etat_zero().

◆ Excision_surf() [2/3]

Lorene::Excision_surf::Excision_surf ( const Excision_surf exc_in)

Copy constructor.

Definition at line 90 of file excision_surf.C.

References set_der_0x0().

◆ Excision_surf() [3/3]

Lorene::Excision_surf::Excision_surf ( FILE *  )

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

◆ ~Excision_surf()

Lorene::Excision_surf::~Excision_surf ( )
virtual

Destructor.

Definition at line 131 of file excision_surf.C.

References del_deriv().

Member Function Documentation

◆ del_deriv()

◆ derive_t_expa()

◆ get_BC_conf_fact_1()

◆ get_BC_conf_fact_2()

◆ get_BC_conf_fact_3()

const Scalar & Lorene::Excision_surf::get_BC_conf_fact_3 ( double  c_theta_lap,
double  c_theta_fin,
Scalar expa_fin 
) const

◆ get_BC_conf_fact_4()

const Scalar & Lorene::Excision_surf::get_BC_conf_fact_4 ( ) const

◆ get_BC_lapse_2()

const Scalar & Lorene::Excision_surf::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.

Definition at line 603 of file excision_surf.C.

References delta_t, Lorene::Scalar::lapang(), lapse, p_get_BC_lapse_2, Lorene::Scalar::set_spectral_va(), Lorene::Scalar::std_spectral_base(), and Lorene::Valeur::ylm().

◆ get_BC_lapse_3()

◆ get_BC_lapse_4()

◆ get_BC_Npsi_2()

const Scalar & Lorene::Excision_surf::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.

Definition at line 1717 of file excision_surf.C.

References conf_fact, delta_t, Lorene::Scalar::lapang(), lapse, p_get_BC_Npsi_2, Lorene::Scalar::set_spectral_va(), Lorene::Scalar::std_spectral_base(), and Lorene::Valeur::ylm().

◆ get_BC_Npsi_3()

const Scalar & Lorene::Excision_surf::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.

Definition at line 1757 of file excision_surf.C.

References conf_fact, Lorene::Map::cost, Lorene::Tensor::get_mp(), lapse, p_get_BC_Npsi_3, Lorene::Scalar::set_spectral_va(), Lorene::sqrt(), and Lorene::Valeur::ylm().

◆ get_BC_Npsi_4()

const Scalar & Lorene::Excision_surf::get_BC_Npsi_4 ( double  Kappa) const

◆ get_BC_Npsi_5()

const Scalar & Lorene::Excision_surf::get_BC_Npsi_5 ( double  Kappa) const

◆ get_BC_shift_2()

◆ get_BC_shift_3()

◆ get_BC_shift_4()

const Vector & Lorene::Excision_surf::get_BC_shift_4 ( Scalar dttheta,
Scalar Ee,
Vector Jj,
Sym_tensor Sij,
double  c_V_lap,
double  epsilon,
bool  sph_sym = true 
) const

◆ get_conf_fact()

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

Returns the conformal factor associated with the surface.

Definition at line 170 of file excision_surf.h.

References conf_fact.

◆ get_delta_t()

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

Returns the timestep used for evolution.

Definition at line 185 of file excision_surf.h.

References delta_t.

◆ get_dt_expa()

const Scalar& Lorene::Excision_surf::get_dt_expa ( ) const
inline

Returns the assumed time derivative of the expansion at t, evolved by functions of this class;.

Definition at line 194 of file excision_surf.h.

References dt_expa.

◆ get_evol_params_from_ID()

void Lorene::Excision_surf::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 a C1 matching with initial data.

Sets also values for expa() and dt_expa() accordingly with initial conditions

Definition at line 195 of file excision_surf.C.

References Lorene::Valeur::c_cf, derive_t_expa(), Lorene::max(), set_dt_expa(), set_expa(), Lorene::Scalar::set_spectral_va(), sph, Lorene::Scalar::std_spectral_base(), Lorene::Spheroid::theta_plus(), Lorene::Mtbl_cf::val_in_bound_jk(), and Lorene::Valeur::ylm().

◆ get_expa()

const Scalar& Lorene::Excision_surf::get_expa ( ) const
inline

Returns the assumed expansion associated to the excised surface at t.

Definition at line 191 of file excision_surf.h.

References expa.

◆ get_gamij()

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

Returns the symmetric tensor $ gamij $.

Definition at line 179 of file excision_surf.h.

References gamij.

◆ get_Kij()

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

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

Definition at line 182 of file excision_surf.h.

References Kij.

◆ get_lapse()

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

Returns the lapse function.

Definition at line 173 of file excision_surf.h.

References lapse.

◆ get_no_of_steps()

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

Returns the internal number of timesteps for one iteration.

Definition at line 188 of file excision_surf.h.

References no_of_steps.

◆ get_shift()

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

Returns the shift vector field.

Definition at line 176 of file excision_surf.h.

References shift.

◆ get_sph()

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

Returns the spheroid.

Definition at line 167 of file excision_surf.h.

References sph.

◆ operator=()

void Lorene::Excision_surf::operator= ( const Excision_surf surf_in)

Assignment to another Excision_surf.

Definition at line 110 of file excision_surf.C.

References conf_fact, del_deriv(), delta_t, dt_expa, expa, gamij, Kij, lapse, no_of_steps, shift, and sph.

◆ sauve()

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

Save in a file.

Definition at line 1851 of file excision_surf.C.

◆ set_conf_fact()

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

Sets the value of the conformal factor.

Definition at line 200 of file excision_surf.h.

References conf_fact, and del_deriv().

◆ set_der_0x0()

◆ set_dt_expa()

Scalar& Lorene::Excision_surf::set_dt_expa ( )
inline

Sets the time derivative of the expansion function on the surface at time t (considering to protect this function)

Definition at line 224 of file excision_surf.h.

References del_deriv(), and dt_expa.

◆ set_expa()

Scalar& Lorene::Excision_surf::set_expa ( )
inline

Sets the expansion function on the surface at time t (considering to protect this function)

Definition at line 220 of file excision_surf.h.

References del_deriv(), and expa.

◆ set_expa_hyperb()

void Lorene::Excision_surf::set_expa_hyperb ( double  alph0,
double  beta0,
double  gamma0 
)

Sets a new value for expansion rescaled over lapse (and its derivative), obtained by hyperbolic evolution.

Parameters for the hyperbolic driver are determined by the function Excision_surf::get_evol_params_from_ID() so that the expansion stays of regularity $C^{1}$ throughout. All manipulated quantities are 2-dimensional.

Definition at line 99 of file set_expa_evol.C.

References Lorene::Scalar::annule_hard(), delta_t, dt_expa, expa, Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Scalar::lapang(), lapse, set_dt_expa(), set_expa(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::val_grid_point(), and Lorene::Valeur::ylm().

◆ set_expa_parab()

void Lorene::Excision_surf::set_expa_parab ( double  c_theta_lap,
double  c_theta_fin,
Scalar expa_fin 
)

◆ set_gamij()

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

Sets the 3d metric of the TimeSlice.

Definition at line 209 of file excision_surf.h.

References del_deriv(), and gamij.

◆ set_Kij()

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

Sets the extrinsic curvature.

Definition at line 212 of file excision_surf.h.

References del_deriv(), and Kij.

◆ set_lapse()

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

Sets the lapse function.

Definition at line 203 of file excision_surf.h.

References del_deriv(), and lapse.

◆ set_shift()

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

Sets the shift vector field.

Definition at line 206 of file excision_surf.h.

References del_deriv(), and shift.

◆ set_sph()

Spheroid& Lorene::Excision_surf::set_sph ( )
inline

Sets a new spheroid from data.

Definition at line 197 of file excision_surf.h.

References del_deriv(), and sph.

Friends And Related Function Documentation

◆ operator<<

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

Display.

Member Data Documentation

◆ conf_fact

Scalar Lorene::Excision_surf::conf_fact
protected

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

Definition at line 53 of file excision_surf.h.

◆ delta_t

double Lorene::Excision_surf::delta_t
protected

The time step for evolution in parabolic drivers.

Definition at line 68 of file excision_surf.h.

◆ dt_expa

Scalar Lorene::Excision_surf::dt_expa
protected

The time derivative of the expansion, derived from Einstein equations and arbitrary evolution.

Definition at line 77 of file excision_surf.h.

◆ expa

Scalar Lorene::Excision_surf::expa
protected

The 2d expansion, directly evolved from the initial excision with Einstein Equations.

Definition at line 74 of file excision_surf.h.

◆ gamij

Metric Lorene::Excision_surf::gamij
protected

The 3-d metric on the slice.

Definition at line 62 of file excision_surf.h.

◆ Kij

Sym_tensor Lorene::Excision_surf::Kij
protected

The 3-d extrinsic curvature on the slice.

Definition at line 65 of file excision_surf.h.

◆ lapse

Scalar Lorene::Excision_surf::lapse
protected

The lapse defined on the 3 slice.

Definition at line 56 of file excision_surf.h.

◆ no_of_steps

double Lorene::Excision_surf::no_of_steps
protected

The internal number of timesteps for one iteration.

Definition at line 71 of file excision_surf.h.

◆ p_derive_t_expa

Scalar* Lorene::Excision_surf::p_derive_t_expa
mutableprotected

Computation of an updated expansion scalar.

Definition at line 94 of file excision_surf.h.

◆ p_get_BC_conf_fact_1

Scalar* Lorene::Excision_surf::p_get_BC_conf_fact_1
mutableprotected

Source of Neumann boundary condition on $ \psi $,.

Definition at line 84 of file excision_surf.h.

◆ p_get_BC_conf_fact_2

Scalar* Lorene::Excision_surf::p_get_BC_conf_fact_2
mutableprotected

Source of Neumann boundary condition on $ \psi $,.

Definition at line 88 of file excision_surf.h.

◆ p_get_BC_conf_fact_3

Scalar* Lorene::Excision_surf::p_get_BC_conf_fact_3
mutableprotected

Source of Neumann boundary condition on $ \psi $,.

Definition at line 89 of file excision_surf.h.

◆ p_get_BC_conf_fact_4

Scalar* Lorene::Excision_surf::p_get_BC_conf_fact_4
mutableprotected

Source of Birichlet boundary condition on $ \psi $,.

Definition at line 90 of file excision_surf.h.

◆ p_get_BC_lapse_1

Scalar* Lorene::Excision_surf::p_get_BC_lapse_1
mutableprotected

Source of Dirichlet boundary condition of $ N $.

Definition at line 85 of file excision_surf.h.

◆ p_get_BC_lapse_2

Scalar* Lorene::Excision_surf::p_get_BC_lapse_2
mutableprotected

Source of Dirichlet boundary condition of $ N $.

Definition at line 91 of file excision_surf.h.

◆ p_get_BC_lapse_3

Scalar* Lorene::Excision_surf::p_get_BC_lapse_3
mutableprotected

Source of Dirichlet condtion on $ N $, based on einstein equations.

Definition at line 92 of file excision_surf.h.

◆ p_get_BC_lapse_4

Scalar* Lorene::Excision_surf::p_get_BC_lapse_4
mutableprotected

Source of Dirichlet condtion on $ N $, based on einstein equations (conservation of isotropic gauge)

Definition at line 93 of file excision_surf.h.

◆ p_get_BC_Npsi_1

Scalar* Lorene::Excision_surf::p_get_BC_Npsi_1
mutableprotected

Source of Neumann boundary condition on $ \psi $.

Definition at line 87 of file excision_surf.h.

◆ p_get_BC_Npsi_2

Scalar* Lorene::Excision_surf::p_get_BC_Npsi_2
mutableprotected

Source of Dirichlet boundary condition on $ N \psi $.

Definition at line 98 of file excision_surf.h.

◆ p_get_BC_Npsi_3

Scalar* Lorene::Excision_surf::p_get_BC_Npsi_3
mutableprotected

Source of Dirichlet boundary condition on $ N \psi $.

Definition at line 99 of file excision_surf.h.

◆ p_get_BC_Npsi_4

Scalar* Lorene::Excision_surf::p_get_BC_Npsi_4
mutableprotected

Source of Dirichlet boundary condition on $ N \psi $.

Definition at line 100 of file excision_surf.h.

◆ p_get_BC_Npsi_5

Scalar* Lorene::Excision_surf::p_get_BC_Npsi_5
mutableprotected

Source of Neumann boundary condition on $ N \psi $.

Definition at line 101 of file excision_surf.h.

◆ p_get_BC_shift_1

Vector* Lorene::Excision_surf::p_get_BC_shift_1
mutableprotected

Source of Dirichlet BC for the shift vector $ \beta^{i} $.

Definition at line 86 of file excision_surf.h.

◆ p_get_BC_shift_2

Vector* Lorene::Excision_surf::p_get_BC_shift_2
mutableprotected

Source of Dirichlet BC for the shift vector $ \beta^{i} $.

Definition at line 95 of file excision_surf.h.

◆ p_get_BC_shift_3

Vector* Lorene::Excision_surf::p_get_BC_shift_3
mutableprotected

Source of Dirichlet BC for the shift vector $ \beta^{i} $, partly derived from kinematical relation.

Definition at line 96 of file excision_surf.h.

◆ p_get_BC_shift_4

Vector* Lorene::Excision_surf::p_get_BC_shift_4
mutableprotected

Source of Dirichlet BC for the shift vector $ \beta^{i} $, partly from projection of Einstein Equations.

Definition at line 97 of file excision_surf.h.

◆ shift

Vector Lorene::Excision_surf::shift
protected

The Shift 3-vector on the slice.

Definition at line 59 of file excision_surf.h.

◆ sph

Spheroid Lorene::Excision_surf::sph
protected

The associated Spheroid object.

Definition at line 50 of file excision_surf.h.


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