LORENE
excision_surf.h
1 /*
2  * Definition of Lorene class Excision_surf, friend class of Spheroid.
3  *
4  */
5 
6 /*
7  * Copyright (c) 2008 Jose-Luis Jaramillo & Nicolas Vasset
8  *
9  * This file is part of LORENE.
10  *
11  * LORENE is free software; you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License version 2
13  * as published by the Free Software Foundation.
14  *
15  * LORENE is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the GNU General Public License
21  * along with LORENE; if not, write to the Free Software
22  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23  *
24  */
25 
26 #ifndef __EXCISIONSURF_H_
27 #define __EXCISIONSURF_H_
28 
29 /*
30  * $Header: /cvsroot/Lorene/C++/Include/excision_surf.h,v 1.7 2014/10/13 08:52:34 j_novak Exp $
31  *
32  */
33 #include "metric.h"
34 #include "spheroid.h"
35 
36 namespace Lorene {
44 
45 
46  // Data :
47  // -----
48  protected:
51 
54 
57 
60 
63 
66 
68  double delta_t;
69 
71  double no_of_steps;
72 
75 
78 
79 
80  // Derived data :
81  // ------------
82  protected:
83 
87  mutable Scalar* p_get_BC_Npsi_1 ;
94  mutable Scalar* p_derive_t_expa ;
98  mutable Scalar* p_get_BC_Npsi_2 ;
99  mutable Scalar* p_get_BC_Npsi_3 ;
102 
103 
104  // Constructors - Destructor
105  // -------------------------
106  public:
107 
119  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) ;
120  Excision_surf(const Excision_surf& ) ;
121 
123  Excision_surf(FILE* ) ;
124 
125  virtual ~Excision_surf() ;
126 
127 
128  // Memory management
129  // -----------------
130  protected:
132  virtual void del_deriv() const ;
133 
135  void set_der_0x0() const ;
136 
137 
138 
139  // Mutators / assignment
140  // ---------------------
141  public:
143  void operator=(const Excision_surf&) ;
144 
145 
146  // Functions/Accessors
147  //---------------------
148 
151  void get_evol_params_from_ID(double alpha, double beta, double gamma, Scalar& Ee, Vector& Jj, Sym_tensor& Ss) ;
152 
154  void set_expa_parab(double c_theta_lap, double c_theta_fin, Scalar& expa_fin) ;
160  void set_expa_hyperb(double alph0, double beta0, double gamma0) ;
161 
162 
163  // Accessors
164  // ---------
165  public:
167  const Spheroid& get_sph() const {return sph; };
168 
170  const Scalar& get_conf_fact() const {return conf_fact; } ;
171 
173  const Scalar& get_lapse() const {return lapse ; } ;
174 
176  const Vector& get_shift() const {return shift ; } ;
177 
179  const Metric& get_gamij() const {return gamij ; } ;
180 
182  const Sym_tensor& get_Kij() const {return Kij ; } ;
183 
185  double get_delta_t() const {return delta_t ;};
186 
188  double get_no_of_steps() const {return no_of_steps ;};
189 
191  const Scalar & get_expa() const {return expa; };
192 
194  const Scalar& get_dt_expa() const {return dt_expa;} ;
195 
197  Spheroid& set_sph() {del_deriv() ; return sph ;};
198 
201 
203  Scalar& set_lapse() {del_deriv() ; return lapse ; } ;
204 
206  Vector& set_shift() {del_deriv() ; return shift ; } ;
207 
209  Metric& set_gamij() {del_deriv() ; return gamij ; } ;
210 
212  Sym_tensor& set_Kij() {del_deriv() ; return Kij ; } ;
213 
214  double set_delta_t() {del_deriv() ; return delta_t ; } ;
215 
216  double set_no_of_steps() {del_deriv() ; return no_of_steps ; } ;
217 
219 
220  Scalar& set_expa() {del_deriv(); return expa;};
221 
223 
225 
226  // Computational functions
227  // -----------------------
228  public:
230  const Scalar& get_BC_conf_fact_1(bool isMOTS = false) const ;
231 // Source for an arbitrary Dirichlet BC on the lapse
232  const Scalar& get_BC_lapse_1(double value) const ;
233 // Source for a global Dirichlet BC on the shift, imposing a conformal Killing symmetry on \f$ \varphi \f$.
234  const Vector& get_BC_shift_1(double Omega) const ;
235 // Source for a Dirichlet arbitrary BC on (N*Psi1)
236  const Scalar& get_BC_Npsi_1(double value) const ;
238  const Scalar& get_BC_conf_fact_2(double c_psi_lap, double c_psi_fin, Scalar& expa_fin) const ;
240  const Scalar& get_BC_conf_fact_3(double c_theta_lap, double c_theta_fin, Scalar& expa_fin) const ;
242  const Scalar& get_BC_conf_fact_4() const ;
244  const Scalar& get_BC_lapse_2(double lapse_fin, double c_lapse_lap, double c_lapse_fi) const ;
246  const Scalar& get_BC_lapse_3(Scalar& dttheta, Scalar& Ee, Vector& Jj, Sym_tensor& Sij, bool sph_sym = true) const ;
248  const Scalar& get_BC_lapse_4 (Scalar& old_nn, Vector& beta_point, Sym_tensor& strain_tens) const ;
250  const Scalar& derive_t_expa(Scalar& Ee, Vector& Jj, Sym_tensor& Sij) const;
252  const Vector& get_BC_shift_2(double c_bb_lap, double c_bb_fin, double c_V_lap, double epsilon) const ;
254  const Vector& get_BC_shift_3(Scalar& dtpsi, double c_V_lap, double epsilon) const ;
256  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 ;
258  const Scalar& get_BC_Npsi_2(double value, double c_npsi_lap, double c_npsi_fin) const ;
260  const Scalar& get_BC_Npsi_3(double n_0, double beta) const ;
262  const Scalar& get_BC_Npsi_4(double Kappa) const ;
264  const Scalar& get_BC_Npsi_5(double Kappa) const ;
265  // Outputs
266  // -------
267  public:
268  virtual void sauve(FILE *) const ;
269 
271  friend ostream& operator<<(ostream& , const Spheroid& ) ;
272 
273 };
274 
275 ostream& operator<<(ostream& , const Spheroid& ) ;
276 
277 
278 }
279 #endif
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) ...
Definition: excision_surf.h:93
Metric & set_gamij()
Sets the 3d metric of the TimeSlice.
Metric for tensor calculation.
Definition: metric.h:90
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...
Definition: set_expa_evol.C:99
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.
Definition: excision_surf.h:68
Sym_tensor Kij
The 3-d extrinsic curvature on the slice.
Definition: excision_surf.h:65
Scalar * p_get_BC_Npsi_2
Source of Dirichlet boundary condition on .
Definition: excision_surf.h:98
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.
Definition: excision_surf.h:74
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.
Definition: excision_surf.h:97
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 ,.
Definition: excision_surf.h:88
Lorene prototypes.
Definition: app_hor.h:67
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).
Definition: scalar.h:393
void operator=(const Excision_surf &)
Assignment to another Excision_surf.
Spheroid sph
The associated Spheroid object.
Definition: excision_surf.h:50
Tensor field of valence 1.
Definition: vector.h:188
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.
Definition: excision_surf.h:92
Scalar * p_get_BC_Npsi_5
Source of Neumann boundary condition on .
Scalar lapse
The lapse defined on the 3 slice.
Definition: excision_surf.h:56
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.
Definition: excision_surf.C:64
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...
Definition: set_expa_evol.C:20
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Definition: spheroid.h:84
Scalar & set_lapse()
Sets the lapse function.
Scalar * p_get_BC_conf_fact_4
Source of Birichlet boundary condition on ,.
Definition: excision_surf.h:90
Scalar * p_get_BC_conf_fact_3
Source of Neumann boundary condition on ,.
Definition: excision_surf.h:89
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 .
Definition: excision_surf.h:99
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.
Definition: excision_surf.h:94
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.
Definition: excision_surf.h:96
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 .
Definition: excision_surf.h:85
double no_of_steps
The internal number of timesteps for one iteration.
Definition: excision_surf.h:71
Surface where boundary conditions for quantities in the bulk will be calculated It relies on geometri...
Definition: excision_surf.h:43
virtual ~Excision_surf()
Destructor.
Scalar dt_expa
The time derivative of the expansion, derived from Einstein equations and arbitrary evolution...
Definition: excision_surf.h:77
Scalar conf_fact
The value of the conformal factor on the 3-slice.
Definition: excision_surf.h:53
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 .
Definition: excision_surf.h:91
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 .
Definition: excision_surf.h:86
Metric gamij
The 3-d metric on the slice.
Definition: excision_surf.h:62
Vector * p_get_BC_shift_2
Source of Dirichlet BC for the shift vector .
Definition: excision_surf.h:95
Scalar * p_get_BC_Npsi_1
Source of Neumann boundary condition on .
Definition: excision_surf.h:87
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.
Definition: excision_surf.h:59
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
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 ,.
Definition: excision_surf.h:84