24char binhor_equations_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Bin_hor/binhor_equations.C,v 1.21 2014/10/13 08:52:42 j_novak Exp $" ;
117#include "utilitaires.h"
127 assert ((relax >0) && (relax<=1)) ;
129 cout <<
"-----------------------------------------------" << endl ;
130 cout <<
"Resolution LAPSE" << endl ;
156 source_un =
hole1.get_psi4()*(
hole1.nn*( aa_quad_un + 0.3333333333333333 *
172 derive_cov(
hole1.ff), 0)
174 .derive_cov(
hole1.ff), 0, 1 ) ;
178 source_un += tmp_un ;
193 source_deux =
hole2.get_psi4()*(
hole2.nn*( aa_quad_deux + 0.3333333333333333
208 derive_cov(
hole2.ff), 0)
210 .derive_cov(
hole2.ff), 0, 1 ) ;
214 source_deux += tmp_deux ;
216 cout <<
"source lapse" << endl <<
norme(source_un) << endl ;
231 lim_un =
hole1.boundary_nn_Dir(lim_nn) ;
232 lim_deux =
hole2.boundary_nn_Dir(lim_nn) ;
234 n_un_temp = n_un_temp - 1./2. ;
235 n_deux_temp = n_deux_temp - 1./2. ;
237 dirichlet_binaire (source_un, source_deux, lim_un, lim_deux,
238 n_un_temp, n_deux_temp, 0, precision) ;
244 lim_un =
hole1.boundary_nn_Neu(lim_nn) ;
245 lim_deux =
hole2.boundary_nn_Neu(lim_nn) ;
247 neumann_binaire (source_un, source_deux, lim_un, lim_deux,
248 n_un_temp, n_deux_temp, 0, precision) ;
253 cout <<
"Unexpected type of boundary conditions for the lapse!"
255 <<
" bound_nn = " << bound_nn << endl ;
262 n_un_temp = n_un_temp + 1./2. ;
263 n_deux_temp = n_deux_temp + 1./2. ;
271 int nz =
hole1.mp.get_mg()->get_nzone() ;
272 cout <<
"lapse auto" << endl <<
norme (n_un_temp) << endl ;
276 "Relative error in the resolution of the equation for the lapse : "
278 for (
int l=0; l<nz; l++) {
279 cout << tdiff_nn(l) <<
" " ;
286 n_un_temp = relax*n_un_temp + (1-relax)*lapse_un_old ;
287 n_deux_temp = relax*n_deux_temp + (1-relax)*lapse_deux_old ;
289 hole1.n_auto = n_un_temp ;
290 hole2.n_auto = n_deux_temp ;
303 assert ((relax>0) && (relax<=1)) ;
305 cout <<
"-----------------------------------------------" << endl ;
306 cout <<
"Resolution PSI" << endl ;
334 tmp_un = 0.125*
hole1.psi_auto * (
hole1.tgam).ricci_scal()
336 .derive_cov(
hole1.ff), 0, 1 ) ;
340 .derive_cov(
hole1.ff), 0) ;
342 source_un = tmp_un -
hole1.psi*
hole1.get_psi4()* ( 0.125* aa_quad_un
363 tmp_deux = 0.125*
hole2.psi_auto * (
hole2.tgam).ricci_scal()
365 .derive_cov(
hole2.ff), 0, 1 ) ;
369 .derive_cov(
hole2.ff), 0) ;
371 source_deux = tmp_deux -
hole2.psi*
hole2.get_psi4()* ( 0.125* aa_quad_deux
376 cout <<
"source psi" << endl <<
norme(source_un) << endl ;
391 lim_un =
hole1.boundary_psi_app_hor() ;
392 lim_deux =
hole2.boundary_psi_app_hor() ;
394 neumann_binaire (source_un, source_deux, lim_un, lim_deux,
395 psi_un_temp, psi_deux_temp, 0, precision) ;
400 cout <<
"Unexpected type of boundary conditions for psi!"
402 <<
" bound_psi = " << bound_psi << endl ;
409 psi_un_temp = psi_un_temp + 1./2. ;
410 psi_deux_temp = psi_deux_temp + 1./2. ;
418 int nz =
hole1.mp.get_mg()->get_nzone() ;
419 cout <<
"psi auto" << endl <<
norme (psi_un_temp) << endl ;
423 "Relative error in the resolution of the equation for psi : "
425 for (
int l=0; l<nz; l++) {
426 cout << tdiff_psi(l) <<
" " ;
433 psi_un_temp = relax*psi_un_temp + (1-relax)*psi_un_old ;
434 psi_deux_temp = relax*psi_deux_temp + (1-relax)*psi_deux_old ;
436 hole1.psi_auto = psi_un_temp ;
437 hole2.psi_auto = psi_deux_temp ;
442 hole1.set_der_0x0() ;
443 hole2.set_der_0x0() ;
456 cout <<
"------------------------------------------------" << endl ;
457 cout <<
"Resolution shift : Omega = " <<
omega << endl ;
487 tmp_vect_un = 2./3.*
hole1.trK.derive_con(
hole1.tgam)
491 source_un += 2.*
hole1.nn * ( tmp_vect_un
493 hole1.aa_auto, 0, 1) ) ;
497 .derive_cov(
hole1.ff), 1, 2)
499 .divergence(
hole1.ff).derive_cov(
hole1.ff), 0)
504 source_un -= vtmp_un ;
506 source_un += 2./3.*
hole1.beta_auto.divergence(
hole1.ff)
530 tmp_vect_deux = 2./3.*
hole2.trK.derive_con(
hole2.tgam)
534 source_deux += 2.*
hole2.nn * ( tmp_vect_deux
540 .derive_cov(
hole2.ff), 1, 2)
542 .divergence(
hole2.ff).derive_cov(
hole2.ff), 0)
547 source_deux -= vtmp_deux ;
549 source_deux += 2./3.*
hole2.beta_auto.divergence(
hole2.ff)
555 Vector source_1 (source_un) ;
556 Vector source_2 (source_deux) ;
560 cout <<
"source shift_x" << endl <<
norme(source_1(1)) << endl ;
561 cout <<
"source shift_y" << endl <<
norme(source_1(2)) << endl ;
562 cout <<
"source shift_z" << endl <<
norme(source_1(3)) << endl ;
565 for (
int i=1 ; i<=3 ; i++) {
577 Valeur lim_x_deux (
hole2.mp.get_mg()-> get_angu()) ;
578 Valeur lim_y_deux (
hole2.mp.get_mg()-> get_angu()) ;
579 Valeur lim_z_deux (
hole2.mp.get_mg()-> get_angu()) ;
581 switch (bound_beta) {
585 lim_x_un =
hole1.boundary_beta_x(
omega, omega_eff) ;
586 lim_y_un =
hole1.boundary_beta_y(
omega, omega_eff) ;
587 lim_z_un =
hole1.boundary_beta_z() ;
589 lim_x_deux =
hole2.boundary_beta_x(
omega, omega_eff) ;
590 lim_y_deux =
hole2.boundary_beta_y(
omega, omega_eff) ;
591 lim_z_deux =
hole2.boundary_beta_z() ;
596 cout <<
"Unexpected type of boundary conditions for beta!"
598 <<
" bound_beta = " << bound_beta << endl ;
614 poisson_vect_binaire (1./3., source_un, source_deux,
615 lim_x_un, lim_y_un, lim_z_un,
616 lim_x_deux, lim_y_deux, lim_z_deux,
617 beta1, beta2, 0, precision) ;
623 for (
int i=1 ; i<=3 ; i++) {
628 cout <<
"shift_auto x" << endl <<
norme(beta1(1)) << endl ;
629 cout <<
"shift_auto y" << endl <<
norme(beta1(2)) << endl ;
630 cout <<
"shift_auto z" << endl <<
norme(beta1(3)) << endl ;
638 int nz =
hole1.mp.get_mg()->get_nzone() ;
643 Tbl tdiff_beta_r =
diffrel(lap_beta(1), source_un(1)) ;
644 Tbl tdiff_beta_t =
diffrel(lap_beta(2), source_un(2)) ;
645 Tbl tdiff_beta_p =
diffrel(lap_beta(3), source_un(3)) ;
648 "Relative error in the resolution of the equation for beta : "
650 cout <<
"r component : " ;
651 for (
int l=0; l<nz; l++) {
652 cout << tdiff_beta_r(l) <<
" " ;
655 cout <<
"t component : " ;
656 for (
int l=0; l<nz; l++) {
657 cout << tdiff_beta_t(l) <<
" " ;
660 cout <<
"p component : " ;
661 for (
int l=0; l<nz; l++) {
662 cout << tdiff_beta_p(l) <<
" " ;
683 if (fabs(
hole1.mp.get_rot_phi()) < 1e-10){
711 if (fabs(
hole2.mp.get_rot_phi()) < 1e-10){
733 beta1_new = relax*(beta1+
hole1.decouple*omdsdp1) + (1-relax)*beta_un_old ;
734 beta2_new = relax*(beta2+
hole2.decouple*omdsdp2) + (1-relax)*beta_deux_old ;
736 hole1.beta_auto = beta1_new ;
737 hole2.beta_auto = beta2_new ;
745 int nnt =
hole1.mp.get_mg()->get_nt(1) ;
746 int nnp =
hole1.mp.get_mg()->get_np(1) ;
750 for (
int k=0; k<nnp; k++)
751 for (
int j=0; j<nnt; j++){
752 if (fabs((
hole1.n_auto+
hole1.n_comp).val_grid_point(1, k, j , 0)) < 1e-8){
759 double diff_un =
hole1.regularisation (
hole1.beta_auto,
761 double diff_deux =
hole2.regularisation (
hole2.beta_auto,
763 hole1.regul = diff_un ;
764 hole2.regul = diff_deux ;
void solve_psi(double precis, double relax, int bound_psi)
Solves the equation for the conformal factor : The fields are the total values excpet those with subs...
double omega
Angular velocity.
Single_hor hole1
Black hole one.
Single_hor hole2
Black hole two.
void solve_shift(double precis, double relax, int bound_beta, double omega_eff)
Solves the equation for the shift, using the Oohara-Nakarmure scheme : The fields are the total value...
void solve_lapse(double precis, double relax, int bound_nn, double lim_nn)
Solves the equation for the lapse : The fields are the total values except those with subscript ,...
Tensor field of valence 0 (or component of a tensorial field).
void raccord(int n)
Performs the matching of the nucleus with respect to the first shell.
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Valeur & set_spectral_va()
Returns va (read/write version).
void annule_hard()
Sets the Scalar to zero in a hard way.
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Class intended to describe valence-2 symmetric tensors.
Values and coefficients of a (real-value) function.
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base ).
Tensor field of valence 1.
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Scalar & set(int)
Read/write access to a component.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
Tensor derive_lie(const Vector &v) const
Computes the Lie derivative of this with respect to some vector field v.
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .