23char bhole_pseudo_kerr_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Bhole/bhole_pseudo_kerr.C,v 1.7 2014/10/13 08:52:40 j_novak Exp $" ;
99 assert ((relax>0) && (relax<=1)) ;
101 cout <<
"Resolution LAPSE" << endl ;
110 for (
int i=0 ; i<3 ; i++) {
111 work.
set() = auxi(i, i) ;
122 Valeur limite (
mp.get_mg()->get_angu()) ;
127 soluce = soluce + 1 ;
130 n_auto.set() = relax*soluce + (1-relax)*lapse_old ;
137 assert ((relax>0) && (relax<=1)) ;
139 cout <<
"Resolution PSI" << endl ;
147 for (
int i=0 ; i<3 ; i++) {
148 work.
set() = auxi(i, i) ;
157 int np =
mp.get_mg()->get_np(1) ;
158 int nt =
mp.get_mg()->get_nt(1) ;
159 Valeur limite (
mp.get_mg()->get_angu()) ;
161 for (
int k=0 ; k<np ; k++)
162 for (
int j=0 ; j<nt ; j++)
167 soluce = soluce + 1 ;
170 psi_auto.set() = relax*soluce + (1-relax)*psi_old ;
178 assert (precision > 0) ;
179 assert ((relax>0) && (relax<=1)) ;
181 cout <<
"resolution SHIFT" << endl ;
192 for (
int i=0 ; i<3 ; i++)
193 if (source(i).get_etat() == ETATQCQ)
201 for (
int i=0 ; i<3 ; i++)
206 int np =
mp.get_mg()->get_np(1) ;
207 int nt =
mp.get_mg()->get_nt(1) ;
209 Mtbl x_mtbl (
mp.get_mg()) ;
211 Mtbl y_mtbl (
mp.get_mg()) ;
217 Base_val** bases =
mp.get_mg()->std_base_vect_cart() ;
219 Valeur lim_x (
mp.get_mg()->get_angu()) ;
221 for (
int k=0 ; k<np ; k++)
222 for (
int j=0 ; j<nt ; j++)
224 lim_x.
base = *bases[0] ;
226 Valeur lim_y (
mp.get_mg()->get_angu()) ;
228 for (
int k=0 ; k<np ; k++)
229 for (
int j=0 ; j<nt ; j++)
230 lim_y.
set(0, k, j, 0) = -
omega*x_mtbl(1, k, j, 0)-
boost[1] ;
231 lim_y.
base = *bases[1] ;
233 Valeur lim_z (
mp.get_mg()->get_angu()) ;
235 for (
int k=0 ; k<np ; k++)
236 for (
int j=0 ; j<nt ; j++)
238 lim_z.
base = *bases[2] ;
241 for (
int i=0 ; i<3 ; i++)
246 poisson_vect_frontiere(1./3., source,
shift_auto, lim_x, lim_y,
247 lim_z, 0, precision, 20) ;
258 for (
int i=0 ; i<3 ; i++) {
263 for (
int i=0 ; i<3 ; i++)
275 Tenseur derive_r (
mp, 1, CON,
mp.get_bvect_cart()) ;
277 for (
int i=0 ; i<3 ; i++)
278 derive_r.
set(i) = tbi(i).dsdr() ;
281 Valeur fonction_radiale (
mp.get_mg()) ;
285 int nz =
mp.get_mg()->get_nzone() ;
286 int np =
mp.get_mg()->get_np(1) ;
287 int nt =
mp.get_mg()->get_nt(1) ;
288 int nr =
mp.get_mg()->get_nr(1) ;
290 double r_0 =
mp.val_r(1, -1, 0, 0) ;
291 double r_1 =
mp.val_r(1, 1, 0, 0) ;
293 for (
int comp=0 ; comp<3 ; comp++) {
295 for (
int k=0 ; k<np ; k++)
296 for (
int j=0 ; j<nt ; j++)
297 for (
int i=0 ; i<nr ; i++)
298 val_hor.
set(1, k, j, i) = derive_r(comp)(1, k, j, 0) ;
300 fonction_radiale =
pow(r_1-
mp.r, 3.)* (
mp.r-r_0)/
pow(r_1-r_0, 3.) ;
301 fonction_radiale.
annule(0) ;
302 fonction_radiale.
annule(2, nz-1) ;
304 enleve = fonction_radiale*val_hor ;
315 if (norm(1) > 1e-5) {
330 Cmp lapse_non_sing (division_xpun(
n_auto(), 0)) ;
334 for (
int i=0 ; i<3 ; i++)
335 for (
int j=i ; j<3 ; j++) {
337 auxi = division_xpun (auxi, 0) ;
338 tkij_auto.set(i, j) = auxi/lapse_non_sing/2. ;
347 double masse =
mp.integrale_surface_infini (integrant) ;
355 double masse =
mp.integrale_surface_infini (integrant) ;
367 for (
int i=0 ; i<3 ; i++)
371 cout <<
"Calcul du moment non valable pour un boost != 0" << endl ;
378 Tenseur vecteur_un (
mp, 1, CON,
mp.get_bvect_cart()) ;
380 for (
int i=0 ; i<3 ; i++)
383 Cmp integrant_un (vecteur_un(0)) ;
386 Tenseur vecteur_deux (
mp, 1, CON,
mp.get_bvect_cart()) ;
388 for (
int i=0 ; i<3 ; i++)
391 Cmp integrant_deux (vecteur_deux(0)) ;
401 double moment =
mp.integrale_surface_infini (-integrant_un+integrant_deux) ;
414 for (
int i=0 ; i<3 ; i++)
418 cout <<
"Calcul du moment non valable pour un boost != 0" << endl ;
434 Tenseur vecteur (
mp, 1, CON,
mp.get_bvect_cart()) ;
436 for (
int i=0 ; i<3 ; i++)
439 vecteur.
annule(
mp.get_mg()->get_nzone()-1) ;
444 double moment =
mp.integrale_surface (integrant,
rayon)/8/M_PI ;
Bases of the spectral expansions.
Tenseur psi_auto
Part of generated by the hole.
Tenseur shift_auto
Part of generated by the hole.
Tenseur_sym tkij_auto
Auto .
double omega
Angular velocity in LORENE's units.
double regul
Intensity of the correction on the shift vector.
Tenseur_sym taij_auto
Part of generated by the hole.
void solve_shift_seul(double precis, double relax)
Solves the equation for ~:
void regularise_seul()
Corrects the shift in the innermost shell, so that it remains and that equals zero on the horizon.
double masse_adm_seul() const
Calculates the ADM mass of the black hole using : .
Tenseur n_auto
Part of N generated by the hole.
double * boost
Lapse on the horizon.
void solve_psi_seul(double relax)
Solves the equation for ~:
double moment_seul_hor() const
Calculates the angular momentum of the black hole using the formula on the horizon : where and H de...
double rayon
Radius of the horizon in LORENE's units.
Map_af & mp
Affine mapping.
void fait_tkij()
Calculates the total .
double masse_komar_seul() const
Calculates the Komar mass of the black hole using : .
void solve_lapse_seul(double relax)
Solves the equation for N ~:
double moment_seul_inf() const
Calculates the angular momentum of the black hole using the formula at infinity : where .
void fait_taij_auto()
Calculates the part of generated by shift_auto .
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
void mult_r_zec()
Multiplication by r in the external compactified domain (ZEC).
Valeur va
The numerical value of the Cmp.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
void annule(int l)
Sets the Cmp to zero in a given domain.
Cmp poisson_dirichlet(const Valeur &limite, int num) const
Is identicall to Cmp::poisson() .
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Cmp poisson_neumann(const Valeur &, int) const
Idem as Cmp::poisson_dirichlet , the boundary condition being on the radial derivative of the solutio...
void raccord(int n)
Performs the matching of the nucleus with respect to the first shell.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void annule(int l)
Sets the Tenseur to zero in a given domain.
void set_std_base()
Set the standard spectal basis of decomposition for each component.
void set_etat_zero()
Sets the logical state to ETATZERO (zero state).
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
int get_etat() const
Returns the logical state.
Values and coefficients of a (real-value) function.
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c ).
void annule(int l)
Sets the Valeur to zero in a given domain.
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
void coef_i() const
Computes the physical value of *this.
const Valeur & mult_st() const
Returns applied to *this.
Base_val base
Bases on which the spectral expansion is performed.
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
void annule_hard()
Sets the Valeur to zero in a hard way.
const Valeur & mult_cp() const
Returns applied to *this.
const Valeur & mult_sp() const
Returns applied to *this.
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Cmp pow(const Cmp &, int)
Power .
Tbl diffrelmax(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (max version).
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
virtual void dsdr(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
Coord xa
Absolute x coordinate.
Coord ya
Absolute y coordinate.