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
|
/*
* Methods for class Et_rot_diff.
*
* (see file et_rot_diff.h for documentation).
*
*/
/*
* Copyright (c) 2001 Eric Gourgoulhon
*
* 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 as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* 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 et_rot_diff_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_diff.C,v 1.4 2014/10/13 08:52:57 j_novak Exp $" ;
/*
* $Id: et_rot_diff.C,v 1.4 2014/10/13 08:52:57 j_novak Exp $
* $Log: et_rot_diff.C,v $
* Revision 1.4 2014/10/13 08:52:57 j_novak
* Lorene classes and functions now belong to the namespace Lorene.
*
* Revision 1.3 2004/03/25 10:29:05 j_novak
* All LORENE's units are now defined in the namespace Unites (in file unites.h).
*
* Revision 1.2 2001/12/04 21:27:53 e_gourgoulhon
*
* All writing/reading to a binary file are now performed according to
* the big endian convention, whatever the system is big endian or
* small endian, thanks to the functions fwrite_be and fread_be
*
* Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
* LORENE
*
* Revision 1.3 2001/10/25 09:20:54 eric
* Ajout de la fonction virtuelle display_poly.
* Affichage de Max nbar, ener et press.
*
* Revision 1.2 2001/10/24 16:23:01 eric
* *** empty log message ***
*
* Revision 1.1 2001/10/19 08:18:10 eric
* Initial revision
*
*
* $Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_diff.C,v 1.4 2014/10/13 08:52:57 j_novak Exp $
*
*/
// Headers C
#include "math.h"
// Headers Lorene
#include "et_rot_diff.h"
#include "eos.h"
#include "nbr_spx.h"
#include "utilitaires.h"
#include "unites.h"
//--------------//
// Constructors //
//--------------//
// Standard constructor
// --------------------
namespace Lorene {
Et_rot_diff::Et_rot_diff(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i,
double (*frot_i)(double, const Tbl&),
double (*primfrot_i)(double, const Tbl&),
const Tbl& par_frot_i)
: Etoile_rot(mp_i, nzet_i, relat, eos_i),
frot(frot_i),
primfrot(primfrot_i),
par_frot(par_frot_i),
omega_field(mp_i),
prim_field(mp_i)
{
// To make sure that omega is not used
omega = __infinity ;
// Initialization to a static state :
omega_field = 0 ;
prim_field = 0 ;
omega_min = 0 ;
omega_max = 0 ;
}
// Copy constructor
// ----------------
Et_rot_diff::Et_rot_diff(const Et_rot_diff& et)
: Etoile_rot(et),
frot(et.frot),
primfrot(et.primfrot),
par_frot(et.par_frot),
omega_field(et.omega_field),
omega_min(et.omega_min),
omega_max(et.omega_max),
prim_field(et.prim_field)
{}
// Constructor from a file
// -----------------------
Et_rot_diff::Et_rot_diff(Map& mp_i, const Eos& eos_i, FILE* fich,
double (*frot_i)(double, const Tbl&),
double (*primfrot_i)(double, const Tbl&) )
: Etoile_rot(mp_i, eos_i, fich),
frot(frot_i),
primfrot(primfrot_i),
par_frot(fich),
omega_field(mp_i),
prim_field(mp_i)
{
Tenseur omega_field_file(mp, fich) ;
omega_field = omega_field_file ;
fait_prim_field() ;
// omega_min and omega_max are read in the file:
fread_be(&omega_min, sizeof(double), 1, fich) ;
fread_be(&omega_max, sizeof(double), 1, fich) ;
}
//------------//
// Destructor //
//------------//
Et_rot_diff::~Et_rot_diff(){}
//--------------//
// Assignment //
//--------------//
// Assignment to another Et_rot_diff
// ---------------------------------
void Et_rot_diff::operator=(const Et_rot_diff& et) {
// Assignment of quantities common to all the derived classes of Etoile_rot
Etoile_rot::operator=(et) ;
// Assignment of proper quantities of class Etoile_rot
frot = et.frot ;
primfrot = et.primfrot ;
par_frot = et.par_frot ;
omega_field = et.omega_field ;
prim_field = et.prim_field ;
omega_min = et.omega_min ;
omega_max = et.omega_max ;
}
//--------------//
// Outputs //
//--------------//
// Save in a file
// --------------
void Et_rot_diff::sauve(FILE* fich) const {
Etoile_rot::sauve(fich) ;
par_frot.sauve(fich) ;
omega_field.sauve(fich) ;
fwrite_be(&omega_min, sizeof(double), 1, fich) ;
fwrite_be(&omega_max, sizeof(double), 1, fich) ;
}
// Printing
// --------
ostream& Et_rot_diff::operator>>(ostream& ost) const {
using namespace Unites ;
Etoile_rot::operator>>(ost) ;
ost << endl ;
ost << "Differentially rotating star" << endl ;
ost << "-----------------------------" << endl ;
ost << endl << "Parameters of F(Omega) : " << endl ;
ost << par_frot << endl ;
ost << "Min, Max of Omega/(2pi) : " << omega_min / (2*M_PI) * f_unit
<< " Hz, " << omega_max / (2*M_PI) * f_unit << " Hz" << endl ;
int lsurf = nzet - 1;
int nt = mp.get_mg()->get_nt(lsurf) ;
int nr = mp.get_mg()->get_nr(lsurf) ;
ost << "Central, equatorial value of Omega: "
<< omega_field()(0, 0, 0, 0) * f_unit << " rad/s, "
<< omega_field()(nzet-1, 0, nt-1, nr-1) * f_unit << " rad/s" << endl ;
ost << "Central, equatorial value of Omega/(2 Pi): "
<< omega_field()(0, 0, 0, 0) * f_unit / (2*M_PI) << " Hz, "
<< omega_field()(nzet-1, 0, nt-1, nr-1) * f_unit / (2*M_PI)
<< " Hz" << endl ;
double nbar_max = max( max( nbar() ) ) ;
double ener_max = max( max( ener() ) ) ;
double press_max = max( max( press() ) ) ;
ost << "Max prop. bar. dens. : " << nbar_max
<< " x 0.1 fm^-3 = " << nbar_max / nbar()(0, 0, 0, 0) << " central"
<< endl ;
ost << "Max prop. ener. dens. (e_max) : " << ener_max
<< " rho_nuc c^2 = " << ener_max / ener()(0, 0, 0, 0) << " central"
<< endl ;
ost << "Max pressure : " << press_max
<< " rho_nuc c^2 = " << press_max / press()(0, 0, 0, 0) << " central"
<< endl ;
// Length scale set by the maximum energy density:
double len_rho = 1. / sqrt( ggrav * ener_max ) ;
ost << endl << "Value of A = par_frot(1) in units of e_max^{1/2} : " <<
par_frot(1) / len_rho << endl ;
ost << "Value of A / r_eq : " <<
par_frot(1) / ray_eq() << endl ;
ost << "r_p/r_eq : " << aplat() << endl ;
ost << "KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
angu_mom() * angu_mom() / pow(len_rho, 0.6666666666666666)
/ pow(mass_b(), 3.3333333333333333)
/ pow(ggrav, 1.3333333333333333) << endl ;
ost << "M e_max^{1/2} : " << ggrav * mass_g() / len_rho << endl ;
ost << "r_eq e_max^{1/2} : " << ray_eq() / len_rho << endl ;
ost << "T/W : " << tsw() << endl ;
ost << "Omega_c / e_max^{1/2} : " << get_omega_c() * len_rho << endl ;
display_poly(ost) ;
return ost ;
}
// display_poly
// ------------
void Et_rot_diff::display_poly(ostream& ost) const {
using namespace Unites ;
Etoile_rot::display_poly( ost ) ;
const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
if (p_eos_poly != 0x0) {
double kappa = p_eos_poly->get_kap() ;
double gamma = p_eos_poly->get_gam() ; ;
// kappa^{n/2}
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 density in terms of rho_unit :
double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
ost.precision(10) ;
ost << " n_max : " << max( max( nbar() ) ) / rho_poly << endl ;
ost << " e_max : " << max( max( ener() ) ) / rho_poly << endl ;
ost << " P_min : " << 2.*M_PI / omega_max / t_poly << endl ;
ost << " P_max : " << 2.*M_PI / omega_min / t_poly << endl ;
int lsurf = nzet - 1;
int nt = mp.get_mg()->get_nt(lsurf) ;
int nr = mp.get_mg()->get_nr(lsurf) ;
ost << " P_eq : " << 2.*M_PI /
omega_field()(nzet-1, 0, nt-1, nr-1) / t_poly << endl ;
}
}
//-----------------------//
// Miscellaneous //
//-----------------------//
double Et_rot_diff::funct_omega(double omeg) const {
return frot(omeg, par_frot) ;
}
double Et_rot_diff::prim_funct_omega(double omeg) const {
return primfrot(omeg, par_frot) ;
}
double Et_rot_diff::get_omega_c() const {
return omega_field()(0, 0, 0, 0) ;
}
}
|