1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443
|
/*
* Methods of class Binary_xcts
* (see file binary_xcts.h for documentation)
*/
/*
* Copyright (c) 2010 Michal Bejger
*
* This file is part of LORENE.
*
* LORENE is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License version 2
* as published by the Free Software Foundation.
*
* LORENE is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with LORENE; if not, write to the Free Software
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*/
char binary_xcts_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binary_xcts/binary_xcts.C,v 1.6 2014/10/13 08:52:45 j_novak Exp $" ;
/*
* $Id: binary_xcts.C,v 1.6 2014/10/13 08:52:45 j_novak Exp $
* $Log: binary_xcts.C,v $
* Revision 1.6 2014/10/13 08:52:45 j_novak
* Lorene classes and functions now belong to the namespace Lorene.
*
* Revision 1.5 2014/10/06 15:12:59 j_novak
* Modified #include directives to use c++ syntax.
*
* Revision 1.4 2010/12/20 09:56:02 m_bejger
* Pointer to the linear momentum added
*
* Revision 1.3 2010/12/09 10:38:46 m_bejger
* virial_vol() added, fait_decouple() removed
*
* Revision 1.2 2010/10/28 12:00:07 m_bejger
* Mass-shedding indicators added to the output in Binary_xcts::write_global
*
* Revision 1.1 2010/05/04 07:35:54 m_bejger
* Initial version
*
* $Header: /cvsroot/Lorene/C++/Source/Binary_xcts/binary_xcts.C,v 1.6 2014/10/13 08:52:45 j_novak Exp $
*
*/
// Headers C
#include <cmath>
// Headers Lorene
#include "binary_xcts.h"
#include "eos.h"
#include "utilitaires.h"
#include "graphique.h"
#include "param.h"
#include "unites.h"
//--------------//
// Constructors //
//--------------//
// Standard constructor
// --------------------
namespace Lorene {
Binary_xcts::Binary_xcts(Map& mp1,
int nzet1,
const Eos& eos1,
int irrot1,
Map& mp2,
int nzet2,
const Eos& eos2,
int irrot2)
: star1(mp1, nzet1, eos1, irrot1),
star2(mp2, nzet2, eos2, irrot2)
{
et[0] = &star1 ;
et[1] = &star2 ;
omega = 0 ;
x_axe = 0 ;
// Pointers of derived quantities initialized to zero :
set_der_0x0() ;
}
// Copy constructor
// ----------------
Binary_xcts::Binary_xcts(const Binary_xcts& bibi)
: star1(bibi.star1),
star2(bibi.star2),
omega(bibi.omega),
x_axe(bibi.x_axe)
{
et[0] = &star1 ;
et[1] = &star2 ;
// Pointers of derived quantities initialized to zero :
set_der_0x0() ;
}
// Constructor from a file
// -----------------------
Binary_xcts::Binary_xcts(Map& mp1,
const Eos& eos1,
Map& mp2,
const Eos& eos2,
FILE* fich)
: star1(mp1, eos1, fich),
star2(mp2, eos2, fich)
{
et[0] = &star1 ;
et[1] = &star2 ;
// omega and x_axe are read in the file:
fread_be(&omega, sizeof(double), 1, fich) ;
fread_be(&x_axe, sizeof(double), 1, fich) ;
// Pointers of derived quantities initialized to zero :
set_der_0x0() ;
}
//------------//
// Destructor //
//------- -----//
Binary_xcts::~Binary_xcts(){
del_deriv() ;
}
//----------------------------------//
// Management of derived quantities //
//----------------------------------//
void Binary_xcts::del_deriv() const {
if (p_mass_adm != 0x0) delete p_mass_adm ;
if (p_mass_kom != 0x0) delete p_mass_kom ;
if (p_angu_mom != 0x0) delete p_angu_mom ;
if (p_lin_mom != 0x0) delete p_lin_mom ;
if (p_total_ener != 0x0) delete p_total_ener ;
if (p_virial != 0x0) delete p_virial ;
if (p_virial_vol != 0x0) delete p_virial_vol ;
set_der_0x0() ;
}
void Binary_xcts::set_der_0x0() const {
p_mass_adm = 0x0 ;
p_mass_kom = 0x0 ;
p_angu_mom = 0x0 ;
p_lin_mom = 0x0 ;
p_total_ener = 0x0 ;
p_virial = 0x0 ;
p_virial_vol = 0x0 ;
}
//--------------//
// Assignment //
//--------------//
// Assignment to another Binary_xcts
// --------------------------------
void Binary_xcts::operator=(const Binary_xcts& bibi) {
star1 = bibi.star1 ;
star2 = bibi.star2 ;
omega = bibi.omega ;
x_axe = bibi.x_axe ;
del_deriv() ; // Deletes all derived quantities
}
//--------------//
// Outputs //
//--------------//
// Save in a file
// --------------
void Binary_xcts::sauve(FILE* fich) const {
star1.sauve(fich) ;
star2.sauve(fich) ;
fwrite_be(&omega, sizeof(double), 1, fich) ;
fwrite_be(&x_axe, sizeof(double), 1, fich) ;
}
// Printing
// --------
ostream& operator<<(ostream& ost, const Binary_xcts& bibi) {
bibi >> ost ;
return ost ;
}
ostream& Binary_xcts::operator>>(ostream& ost) const {
using namespace Unites ;
ost << endl ;
ost << "Binary neutron stars" << endl ;
ost << "=============" << endl ;
ost << endl <<
"Orbital angular velocity : " << omega * f_unit << " rad/s" << endl ;
ost << endl <<
"Coordinate separation between the two stellar centers : "
<< separation() / km << " km" << endl ;
ost <<
"Absolute coordinate X of the rotation axis : " << x_axe / km
<< " km" << endl ;
ost << endl << "Star 1 : " << endl ;
ost << "====== " << endl ;
ost << star1 << endl ;
ost << "Star 2 : " << endl ;
ost << "====== " << endl ;
ost << star2 << endl ;
return ost ;
}
// Display in polytropic units
// ---------------------------
void Binary_xcts::display_poly(ostream& ost) const {
using namespace Unites ;
const Eos* p_eos1 = &( star1.get_eos() ) ;
const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;
if (p_eos_poly != 0x0) {
assert( star1.get_eos() == star2.get_eos() ) ;
double kappa = p_eos_poly->get_kap() ;
double gamma = p_eos_poly->get_gam() ; ;
double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
// Polytropic unit of length in terms of r_unit :
double r_poly = kap_ns2 / sqrt(ggrav) ;
// Polytropic unit of time in terms of t_unit :
double t_poly = r_poly ;
// Polytropic unit of mass in terms of m_unit :
double m_poly = r_poly / ggrav ;
// Polytropic unit of angular momentum in terms of j_unit :
double j_poly = r_poly * r_poly / ggrav ;
ost.precision(10) ;
ost << endl << "Quantities in polytropic units : " << endl ;
ost << "==============================" << endl ;
ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
ost << " d_e_max : " << separation() / r_poly << endl ;
ost << " d_G : "
<< ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly
<< endl ;
ost << " Omega : " << omega * t_poly << endl ;
ost << " J : " << angu_mom()(2) / j_poly << endl ;
ost << " M_ADM : " << mass_adm() / m_poly << endl ;
ost << " M_Komar : " << mass_kom() / m_poly << endl ;
ost << " M_bar(star 1) : " << star1.mass_b() / m_poly << endl ;
ost << " M_bar(star 2) : " << star2.mass_b() / m_poly << endl ;
ost << " R_0(star 1) : " <<
0.5 * ( star1.ray_eq() + star1.ray_eq_pi() ) / r_poly << endl ;
ost << " R_0(star 2) : " <<
0.5 * ( star2.ray_eq() + star2.ray_eq_pi() ) / r_poly << endl ;
}
}
void Binary_xcts::write_global(ostream& ost) const {
using namespace Unites ;
const Map& mp1 = star1.get_mp() ;
const Mg3d* mg1 = mp1.get_mg() ;
int nz1 = mg1->get_nzone() ;
// Mass-shedding indicators
double dent1_eq = (star1.ent).dsdr().val_point(star1.ray_eq(),M_PI/2.,0.) ;
double dent1_pole = (star1.ent).dsdr().val_point(star1.ray_pole(),0.,0.) ;
double chi1 = fabs( dent1_eq / dent1_pole ) ;
double dent2_eq = (star2.ent).dsdr().val_point(star2.ray_eq(),M_PI/2.,0.) ;
double dent2_pole = (star2.ent).dsdr().val_point(star2.ray_pole(),0.,0.) ;
double chi2 = fabs( dent2_eq / dent2_pole ) ;
ost.precision(5) ;
ost << "# Grid 1 : " << nz1 << "x"
<< mg1->get_nr(0) << "x" << mg1->get_nt(0) << "x" << mg1->get_np(0)
<< " R_out(l) [km] : " ;
for (int l=0; l<nz1; l++) {
ost << " " << mp1.val_r(l, 1., M_PI/2, 0) / km ;
}
ost << endl ;
ost << "# VE(M) VE(M) [vol]" << endl ;
ost.setf(ios::scientific) ;
ost.width(14) ;
ost << virial() << " " << virial_vol() << endl ;
ost << "# d [km] "
<< " d_G [km] "
<< " d/(a1 +a1') "
<< " f [Hz] "
<< " M_ADM [M_sol] "
<< " M_ADM_vol [M_sol] "
<< " M_Komar [M_sol] "
<< " M_Komar_vol [M_sol] "
<< " J [G M_sol^2/c] " << endl ;
ost.precision(14) ;
ost.width(20) ;
ost << separation() / km ; ost.width(22) ;
ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / km ; ost.width(22) ;
ost << separation() / (star1.ray_eq() + star2.ray_eq()) ; ost.width(22) ;
ost << omega / (2*M_PI)* f_unit ; ost.width(22) ;
ost << mass_adm() / msol ; ost.width(22) ;
ost << mass_adm_vol() / msol ; ost.width(22) ;
ost << mass_kom() / msol ; ost.width(22) ;
ost << mass_kom_vol() / msol ; ost.width(22) ;
ost << angu_mom()(2)/ ( qpig / (4* M_PI) * msol*msol) << endl ;
ost << "# H_c(1)[c^2] "
<< " e_c(1)[rho_nuc] "
<< " M_B(1) [M_sol] "
<< " r_eq(1) [km] "
<< " a2/a1(1) "
<< " a3/a1(1) "
<< " chi1 " << endl ;
ost.width(20) ;
ost << star1.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
ost << star1.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
ost << star1.mass_b() / msol ; ost.width(22) ;
ost << star1.ray_eq() / km ; ost.width(22) ;
ost << star1.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
ost << star1.ray_pole() / star1.ray_eq() ; ost.width(22) ;
ost << chi1 << endl ;
ost << "# H_c(2)[c^2] "
<< " e_c(2)[rho_nuc] "
<< " M_B(2) [M_sol] "
<< " r_eq(2) [km] "
<< " a2/a1(2) "
<< " a3/a1(2) "
<< " chi2 " << endl ;
ost.width(20) ;
ost << star2.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
ost << star2.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
ost << star2.mass_b() / msol ; ost.width(22) ;
ost << star2.ray_eq() / km ; ost.width(22) ;
ost << star2.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
ost << star2.ray_pole() / star1.ray_eq() ; ost.width(22) ;
ost << chi2 << endl ;
// Quantities in polytropic units if the EOS is a polytropic one
// -------------------------------------------------------------
const Eos* p_eos1 = &( star1.get_eos() ) ;
const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;
if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
double kappa = p_eos_poly->get_kap() ;
double gamma = p_eos_poly->get_gam() ; ;
double kap_ns2 = pow( kappa, 0.5 /(gamma-1.) ) ;
// Polytropic unit of length in terms of r_unit :
double r_poly = kap_ns2 / sqrt(ggrav) ;
// Polytropic unit of time in terms of t_unit :
double t_poly = r_poly ;
// Polytropic unit of mass in terms of m_unit :
double m_poly = r_poly / ggrav ;
// Polytropic unit of angular momentum in terms of j_unit :
double j_poly = r_poly * r_poly / ggrav ;
ost << "# d [poly] "
<< " d_G [poly] "
<< " Omega [poly] "
<< " M_ADM [poly] "
<< " J [poly] "
<< " M_B(1) [poly] "
<< " M_B(2) [poly] " << endl ;
ost.width(20) ;
ost << separation() / r_poly ; ost.width(22) ;
ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly ; ost.width(22) ;
ost << omega * t_poly ; ost.width(22) ;
ost << mass_adm() / m_poly ; ost.width(22) ;
ost << angu_mom()(2) / j_poly ; ost.width(22) ;
ost << star1.mass_b() / m_poly ; ost.width(22) ;
ost << star2.mass_b() / m_poly << endl ;
}
}
//-------------------------------//
// Miscellaneous //
//-------------------------------//
double Binary_xcts::separation() const {
double dx = star1.mp.get_ori_x() - star2.mp.get_ori_x() ;
double dy = star1.mp.get_ori_y() - star2.mp.get_ori_y() ;
double dz = star1.mp.get_ori_z() - star2.mp.get_ori_z() ;
return sqrt( dx*dx + dy*dy + dz*dz ) ;
}
}
|