00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 char bin_hor_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bin_hor/bin_hor.C,v 1.10 2007/04/13 15:28:55 f_limousin Exp $" ;
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075 #include <stdlib.h>
00076 #include <math.h>
00077
00078
00079 #include "nbr_spx.h"
00080 #include "tenseur.h"
00081 #include "tensor.h"
00082 #include "isol_hor.h"
00083 #include "proto.h"
00084 #include "utilitaires.h"
00085
00086
00087
00088
00089
00090 Bin_hor::Bin_hor (Map_af& mp1, Map_af& mp2) :
00091 hole1(mp1), hole2(mp2), omega(0){
00092
00093 holes[0] = &hole1 ;
00094 holes[1] = &hole2 ;
00095 }
00096
00097
00098
00099
00100 Bin_hor::Bin_hor (const Bin_hor& source) :
00101 hole1(source.hole1), hole2(source.hole2), omega(source.omega) {
00102
00103 holes[0] = &hole1 ;
00104 holes[1] = &hole2 ;
00105 }
00106
00107
00108
00109
00110 Bin_hor::Bin_hor(Map_af& mp1, Map_af& mp2, FILE* fich)
00111 : hole1(mp1, fich),
00112 hole2(mp2, fich),
00113 omega(0) {
00114
00115 fread_be(&omega, sizeof(double), 1, fich) ;
00116 holes[0] = &hole1 ;
00117 holes[1] = &hole2 ;
00118
00119 }
00120
00121
00122
00123
00124
00125 Bin_hor::~Bin_hor () {
00126 }
00127
00128
00129
00130
00131
00132 void Bin_hor::operator= (const Bin_hor& source) {
00133 hole1 = source.hole1 ;
00134 hole2 = source.hole2 ;
00135
00136 omega = source.omega ;
00137 }
00138
00139
00140
00141
00142
00143
00144 void Bin_hor::sauve(FILE* fich) const {
00145
00146 hole1.sauve(fich) ;
00147 hole2.sauve(fich) ;
00148 fwrite_be (&omega, sizeof(double), 1, fich) ;
00149
00150 }
00151
00152
00153
00154 void Bin_hor::init_bin_hor() {
00155 set_omega (0) ;
00156 hole1.init_bhole() ;
00157 hole2.init_bhole() ;
00158
00159 hole1.psi_comp_import(hole2) ;
00160 hole2.psi_comp_import(hole1) ;
00161
00162 hole1.n_comp_import(hole2) ;
00163 hole2.n_comp_import(hole1) ;
00164
00165 decouple() ;
00166 extrinsic_curvature() ;
00167
00168 }
00169
00170
00171 void Bin_hor::write_global(ostream& ost, double lim_nn, int bound_nn,
00172 int bound_psi, int bound_beta, double alpha) const {
00173
00174 double distance = hole1.get_mp().get_ori_x() - hole2.get_mp().get_ori_x() ;
00175 double mass_adm = adm_mass() ;
00176 double mass_komar = komar_mass() ;
00177 double mass_area = sqrt(hole1.area_hor()/16/M_PI) +
00178 sqrt(hole2.area_hor()/16/M_PI) ;
00179 double J_adm = ang_mom_adm() ;
00180 double J_hor = ang_mom_hor() ;
00181 double j1 = hole1.ang_mom_hor() ;
00182 double j2 = hole2.ang_mom_hor() ;
00183 double mass_ih1 = hole1.mass_hor() ;
00184 double mass_ih2 = hole2.mass_hor() ;
00185 double mass_ih = mass_ih1 + mass_ih2 ;
00186 double omega1 = hole1.omega_hor() ;
00187 double omega2 = hole2.omega_hor() ;
00188
00189
00190
00191
00192
00193 double orientation1 = hole1.mp.get_rot_phi() ;
00194 assert ((orientation1 == 0) || (orientation1 == M_PI)) ;
00195 int aligne1 = (orientation1 == 0) ? 1 : -1 ;
00196
00197 Vector angular1 (hole1.mp, CON, hole1.mp.get_bvect_cart()) ;
00198 Scalar yya1 (hole1.mp) ;
00199 yya1 = hole1.mp.ya ;
00200 Scalar xxa1 (hole1.mp) ;
00201 xxa1 = hole1.mp.xa ;
00202
00203 angular1.set(1) = aligne1 * omega * yya1 ;
00204 angular1.set(2) = - aligne1 * omega * xxa1 ;
00205 angular1.set(3).annule_hard() ;
00206
00207 angular1.set(1).set_spectral_va()
00208 .set_base(*(hole1.mp.get_mg()->std_base_vect_cart()[0])) ;
00209 angular1.set(2).set_spectral_va()
00210 .set_base(*(hole1.mp.get_mg()->std_base_vect_cart()[1])) ;
00211 angular1.set(3).set_spectral_va()
00212 .set_base(*(hole1.mp.get_mg()->std_base_vect_cart()[2])) ;
00213
00214 angular1.change_triad(hole1.mp.get_bvect_spher()) ;
00215
00216
00217 double orientation2 = hole2.mp.get_rot_phi() ;
00218 assert ((orientation2 == 0) || (orientation2 == M_PI)) ;
00219 int aligne2 = (orientation2 == 0) ? 1 : -1 ;
00220
00221 Vector angular2 (hole2.mp, CON, hole2.mp.get_bvect_cart()) ;
00222 Scalar yya2 (hole2.mp) ;
00223 yya2 = hole2.mp.ya ;
00224 Scalar xxa2 (hole2.mp) ;
00225 xxa2 = hole2.mp.xa ;
00226
00227 angular2.set(1) = aligne2 * omega * yya2 ;
00228 angular2.set(2) = - aligne2 * omega * xxa2 ;
00229 angular2.set(3).annule_hard() ;
00230
00231 angular2.set(1).set_spectral_va()
00232 .set_base(*(hole2.mp.get_mg()->std_base_vect_cart()[0])) ;
00233 angular2.set(2).set_spectral_va()
00234 .set_base(*(hole2.mp.get_mg()->std_base_vect_cart()[1])) ;
00235 angular2.set(3).set_spectral_va()
00236 .set_base(*(hole2.mp.get_mg()->std_base_vect_cart()[2])) ;
00237
00238 angular2.change_triad(hole2.mp.get_bvect_spher()) ;
00239
00240
00241 Scalar btilde1 (hole1.b_tilde() -
00242 contract(angular1, 0, hole1.tgam.radial_vect().up_down(hole1.tgam), 0)) ;
00243 Scalar btilde2 (hole2.b_tilde() -
00244 contract(angular2, 0, hole2.tgam.radial_vect().up_down(hole2.tgam), 0)) ;
00245
00246
00247
00248
00249 Vector integrand_un (hole1.mp, COV, hole1.mp.get_bvect_spher()) ;
00250 integrand_un = hole1.nn.derive_cov(hole1.ff)*pow(hole1.psi, 2)
00251 - btilde1*contract(hole1.get_k_dd(), 1,
00252 hole1.tgam.radial_vect(), 0)*pow(hole1.psi, 2) ;
00253 integrand_un.std_spectral_base() ;
00254
00255 Vector integrand_deux (hole2.mp, COV, hole2.mp.get_bvect_spher()) ;
00256 integrand_deux = hole2.nn.derive_cov(hole2.ff)*pow(hole2.psi, 2)
00257 - btilde2*contract(hole2.get_k_dd(), 1,
00258 hole2.tgam.radial_vect(), 0)*pow(hole2.psi, 2) ;
00259 integrand_deux.std_spectral_base() ;
00260
00261 double horizon = hole1.mp.integrale_surface(integrand_un(1),
00262 hole1.get_radius())+
00263 hole2.mp.integrale_surface(integrand_deux(1), hole2.get_radius()) ;
00264
00265 horizon /= 4*M_PI ;
00266
00267 double J_smarr = (mass_komar - horizon) / 2. / omega ;
00268
00269 ost.precision(8) ;
00270 ost << "# Grid : " << hole1.mp.get_mg()->get_nr(1) << "x"
00271 << hole1.mp.get_mg()->get_nt(1) << "x"
00272 << hole1.mp.get_mg()->get_np(1) << " R_out(l) : " ;
00273
00274 for (int ll=0; ll<hole1.mp.get_mg()->get_nzone(); ll++) {
00275 ost << " " << hole1.mp.val_r(ll, 1., M_PI/2, 0) ;
00276 }
00277 ost << endl ;
00278 ost << "# bound N, lim N : " << bound_nn << " " << lim_nn
00279 << " - bound Psi : " << bound_psi << " - bound shift : " << bound_beta
00280 << " alpha = " << alpha << endl ;
00281
00282 ost << "# distance omega Mass_ADM Mass_K M_area J_ADM J_hor" << endl ;
00283 ost << distance << " " ;
00284 ost << omega << " " ;
00285 ost << mass_adm << " " ;
00286 ost << mass_komar << " " ;
00287 ost << mass_area << " " ;
00288 ost << J_adm << " " ;
00289 ost << J_hor << endl ;
00290 ost << "# mass_ih1 mass_ih2 mass_ih j1 J2 omega1 omega2" << endl ;
00291 ost << mass_ih1 << " " ;
00292 ost << mass_ih2 << " " ;
00293 ost << mass_ih << " " ;
00294 ost << j1 << " " ;
00295 ost << j2 << " " ;
00296 ost << omega1 << " " ;
00297 ost << omega2 << endl ;
00298 ost << "# ADM_mass/M_area J/M_area2 omega*M_area" << endl ;
00299 ost << mass_adm / mass_area << " " ;
00300 ost << J_adm /mass_area / mass_area << " " ;
00301 ost << omega * mass_area << endl ;
00302 ost << "# Diff J_hor/J_ADM Diff J_ADM/J_Smarr Diff J_hor/J_smarr"
00303 << endl ;
00304 ost << fabs(J_adm - J_hor) / J_adm << " " << fabs(J_adm - J_smarr) / J_adm
00305 << " " << fabs(J_hor - J_smarr) / J_hor << endl ;
00306
00307
00308 }
00309