LORENE
star_rot.h
1/*
2 * Definition of Lorene class Star_rot
3 *
4 */
5
6/*
7 * Copyright (c) 2010 Eric Gourgoulhon
8 *
9 * This file is part of LORENE.
10 *
11 * LORENE is free software; you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License as published by
13 * the Free Software Foundation; either version 2 of the License, or
14 * (at your option) any later version.
15 *
16 * LORENE is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 * GNU General Public License for more details.
20 *
21 * You should have received a copy of the GNU General Public License
22 * along with LORENE; if not, write to the Free Software
23 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
24 *
25 */
26
27
28#ifndef __STAR_ROT_H_
29#define __STAR_ROT_H_
30
31/*
32 * $Id: star_rot.h,v 1.10 2023/07/04 09:01:09 j_novak Exp $
33 * $Log: star_rot.h,v $
34 * Revision 1.10 2023/07/04 09:01:09 j_novak
35 * Changed the name r_pol_circ -> r_circ_merid to explicitly show the radius is a circumferential one, along a meridian.
36 *
37 * Revision 1.9 2023/07/04 07:28:34 j_novak
38 * Added method "r_pol_circ()", to compute polar circumferential radius, with a circumference along a meridian.
39 *
40 * Revision 1.8 2017/10/20 13:55:41 j_novak
41 * Calss Star_rot is now friend of class Star.
42 *
43 * Revision 1.7 2017/04/11 10:46:55 m_bejger
44 * Star_rot::surf_grav() - surface gravity values along the theta direction
45 *
46 * Revision 1.6 2015/05/19 09:30:55 j_novak
47 * New methods for computing the area of the star and its mean radius.
48 *
49 * Revision 1.5 2014/10/13 08:52:36 j_novak
50 * Lorene classes and functions now belong to the namespace Lorene.
51 *
52 * Revision 1.4 2010/02/08 10:56:30 j_novak
53 * Added a few things missing for the reading from resulting file.
54 *
55 * Revision 1.3 2010/02/02 13:35:00 e_gourgoulhon
56 * Remove the "under construction" mark.
57 *
58 * Revision 1.2 2010/01/25 18:14:05 e_gourgoulhon
59 * Added member unsurc2
60 * Added method is_relativistic()
61 * Suppressed method f_eccentric
62 *
63 * Revision 1.1 2010/01/24 16:07:45 e_gourgoulhon
64 * New class Star_rot.
65 *
66 *
67 * $Header: /cvsroot/Lorene/C++/Include/star_rot.h,v 1.10 2023/07/04 09:01:09 j_novak Exp $
68 *
69 */
70
71// Headers Lorene
72#include "star.h"
73
74namespace Lorene {
75class Eos ;
76
77 //--------------------------//
78 // class Star_rot //
79 //--------------------------//
80
98class Star_rot : public Star {
99
100 // Data :
101 // -----
102 protected:
107
111 double unsurc2 ;
112
113 double omega ;
114
117
120
123
126
131
134
139
144
147
150
163
173
180
199
205
211
216
221
229
238
239
240 // Derived data :
241 // ------------
242 protected:
243
244 mutable double* p_angu_mom ;
245 mutable double* p_tsw ;
246 mutable double* p_grv2 ;
247 mutable double* p_grv3 ;
248 mutable double* p_r_circ ;
249 mutable double* p_r_circ_merid ;
250 mutable double* p_aplat ;
251 mutable double* p_area ;
252 mutable double* p_z_eqf ;
253 mutable double* p_z_eqb ;
254 mutable double* p_z_pole ;
255 mutable double* p_mom_quad ;
256 mutable double* p_r_isco ;
257 mutable double* p_f_isco ;
259 mutable double* p_espec_isco ;
261 mutable double* p_lspec_isco ;
262 mutable double* p_f_eq ;
263 mutable Tbl* p_surf_grav ;
264
265
266
267 // Constructors - Destructor
268 // -------------------------
269 public:
278 Star_rot(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
279
280
281 Star_rot(const Star_rot& ) ;
282
290 Star_rot(Map& mp_i, const Eos& eos_i, FILE* fich) ;
291
292 virtual ~Star_rot() ;
293
294
295 // Memory management
296 // -----------------
297 protected:
299 virtual void del_deriv() const ;
300
302 virtual void set_der_0x0() const ;
303
307 virtual void del_hydro_euler() ;
308
309
310 // Mutators / assignment
311 // ---------------------
312 public:
314 void operator=(const Star_rot& ) ;
315
316 // Accessors
317 // ---------
318 public:
322 bool is_relativistic() const {return relativistic; } ;
323
327 virtual double get_omega_c() const ;
328
330 const Scalar& get_bbb() const {return bbb;} ;
331
333 const Scalar& get_a_car() const {return a_car;} ;
334
336 const Scalar& get_b_car() const {return b_car;} ;
337
339 const Scalar& get_nphi() const {return nphi;} ;
340
344 const Scalar& get_tnphi() const {return tnphi;} ;
345
347 const Scalar& get_uuu() const {return uuu;} ;
348
352 const Scalar& get_nuf() const {return nuf;} ;
353
357 const Scalar& get_nuq() const {return nuq;} ;
358
360 const Scalar& get_dzeta() const {return dzeta;} ;
361
363 const Scalar& get_tggg() const {return tggg;} ;
364
377 const Vector& get_w_shift() const {return w_shift;} ;
378
391 const Scalar& get_khi_shift() const {return khi_shift;} ;
392
398 const Sym_tensor& get_tkij() const {return tkij;} ;
399
417 const Scalar& get_ak_car() const {return ak_car;} ;
418
419 // Outputs
420 // -------
421 public:
422 virtual void sauve(FILE* ) const ;
423
425 virtual void display_poly(ostream& ) const ;
426
427 protected:
429 virtual ostream& operator>>(ostream& ) const ;
430
432 virtual void partial_display(ostream& ) const ;
433
434 // Global quantities
435 // -----------------
436 public:
437
445 virtual const Itbl& l_surf() const ;
446
447 virtual double mass_b() const ;
448 virtual double mass_g() const ;
449 virtual double angu_mom() const ;
450 virtual double tsw() const ;
451
455 virtual double grv2() const ;
456
468 virtual double grv3(ostream* ost = 0x0) const ;
469
470 virtual double r_circ() const ;
471 virtual double r_circ_merid() const ;
472 virtual double aplat() const ;
473 virtual double area() const ;
475 virtual double mean_radius() const ;
476 virtual double z_eqf() const ;
477 virtual double z_eqb() const ;
478 virtual double z_pole() const ;
479
489 virtual double mom_quad() const ;
490
494 virtual const Tbl& surf_grav() const ;
495
502 virtual double r_isco(ostream* ost = 0x0) const ;
503
505 virtual double f_isco() const ;
506
508 virtual double espec_isco() const ;
509
511 virtual double lspec_isco() const ;
512
514 virtual double f_eq() const ;
515
516
517 // Computational routines
518 // ----------------------
519 public:
530 virtual void hydro_euler() ;
531
542 void update_metric() ;
543
552 void fait_shift() ;
553
557 void fait_nphi() ;
558
562 void extrinsic_curvature() ;
563
593 static double lambda_grv2(const Scalar& sou_m, const Scalar& sou_q) ;
594
674 virtual void equilibrium(double ent_c, double omega0, double fact_omega,
675 int nzadapt, const Tbl& ent_limit,
676 const Itbl& icontrol, const Tbl& control,
677 double mbar_wanted, double aexp_mass,
678 Tbl& diff, Param* = 0x0) ;
679
680
681 friend class Star ;
682};
683
684}
685#endif
Equation of state base class.
Definition eos.h:206
Basic integer array class.
Definition itbl.h:122
Parameter storage.
Definition param.h:125
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
double * p_angu_mom
Angular momentum.
Definition star_rot.h:244
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition star_rot.C:675
virtual double lspec_isco() const
Angular momentum of a particle on the ISCO.
virtual double mean_radius() const
Mean star radius from the area .
double * p_r_isco
Circumferential radius of the ISCO.
Definition star_rot.h:256
double * p_aplat
Flatening r_pole/r_eq.
Definition star_rot.h:250
Star_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition star_rot.C:97
double * p_mom_quad
Quadrupole moment.
Definition star_rot.h:255
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Definition star_rot.h:220
Sym_tensor tkij
Tensor related to the extrinsic curvature tensor by .
Definition star_rot.h:179
virtual double espec_isco() const
Energy of a particle on the ISCO.
virtual const Tbl & surf_grav() const
Surface gravity (table along the theta direction).
virtual double mom_quad() const
Quadrupole moment.
double * p_z_eqf
Forward redshift factor at equator.
Definition star_rot.h:252
Scalar tggg
Metric potential .
Definition star_rot.h:149
Scalar tnphi
Component of the shift vector.
Definition star_rot.h:130
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog.
Definition star_rot.C:730
virtual ~Star_rot()
Destructor.
Definition star_rot.C:300
virtual double z_pole() const
Redshift factor at North pole.
const Sym_tensor & get_tkij() const
Returns the tensor related to the extrinsic curvature tensor by .
Definition star_rot.h:398
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
virtual double z_eqb() const
Backward redshift factor at equator.
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition star_rot.h:111
const Scalar & get_dzeta() const
Returns the Metric potential .
Definition star_rot.h:360
const Scalar & get_tggg() const
Returns the Metric potential .
Definition star_rot.h:363
void extrinsic_curvature()
Computes tkij and ak_car from shift , nnn and b_car .
Vector w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition star_rot.h:162
const Vector & get_w_shift() const
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition star_rot.h:377
Scalar b_car
Square of the metric factor B.
Definition star_rot.h:122
void operator=(const Star_rot &)
Assignment to another Star_rot.
Definition star_rot.C:376
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition star_rot.h:261
Scalar bbb
Metric factor B.
Definition star_rot.h:119
virtual double mass_b() const
Baryon mass.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] ).
Definition star_rot.C:666
const Scalar & get_khi_shift() const
Returns the scalar used in the decomposition of shift following Shibata's prescription [Prog.
Definition star_rot.h:391
double omega
Rotation angular velocity ([f_unit] ).
Definition star_rot.h:113
const Scalar & get_nuq() const
Returns the Part of the Metric potential = logn generated by the quadratic terms.
Definition star_rot.h:357
Scalar ak_car
Scalar .
Definition star_rot.h:198
virtual double r_circ() const
Circumferential equatorial radius.
Scalar uuu
Norm of u_euler.
Definition star_rot.h:133
virtual double z_eqf() const
Forward redshift factor at equator.
const Scalar & get_ak_car() const
Returns the scalar .
Definition star_rot.h:417
Scalar nphi
Metric coefficient .
Definition star_rot.h:125
const Scalar & get_nphi() const
Returns the metric coefficient .
Definition star_rot.h:339
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition star_rot.C:337
double * p_f_isco
Orbital frequency of the ISCO.
Definition star_rot.h:257
void update_metric()
Computes metric coefficients from known potentials.
virtual double f_isco() const
Orbital frequency at the innermost stable circular orbit (ISCO).
virtual double mass_g() const
Gravitational mass.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition star_rot.C:445
virtual void sauve(FILE *) const
Save in a file.
Definition star_rot.C:418
double * p_area
Integrated surface area.
Definition star_rot.h:251
virtual double aplat() const
Flatening r_pole/r_eq.
double * p_grv3
Error on the virial identity GRV3.
Definition star_rot.h:247
Scalar ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition star_rot.h:210
const Scalar & get_a_car() const
Returns the square of the metric factor A.
Definition star_rot.h:333
virtual const Itbl & l_surf() const
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on...
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition star_rot.h:259
double * p_tsw
Ratio T/W.
Definition star_rot.h:245
const Scalar & get_b_car() const
Returns the square of the metric factor B.
Definition star_rot.h:336
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
Definition star_rot.h:138
Tbl * p_surf_grav
Surface gravity (along the theta direction).
Definition star_rot.h:263
bool is_relativistic() const
Returns true for a relativistic star, false for a Newtonian one.
Definition star_rot.h:322
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition star_rot.h:106
double * p_f_eq
Orbital frequency at the equator.
Definition star_rot.h:262
Scalar ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta .
Definition star_rot.h:215
virtual double tsw() const
Ratio T/W.
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Definition star_rot.C:362
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition star_rot.h:228
Scalar dzeta
Metric potential .
Definition star_rot.h:146
virtual double area() const
Integrated surface area in .
Scalar a_car
Square of the metric factor A.
Definition star_rot.h:116
Scalar ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition star_rot.h:204
Scalar nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition star_rot.h:143
const Scalar & get_uuu() const
Returns the norm of u_euler.
Definition star_rot.h:347
virtual double r_circ_merid() const
Crcumferential meridional radius.
double * p_grv2
Error on the virial identity GRV2.
Definition star_rot.h:246
virtual void del_deriv() const
Deletes all the derived quantities.
Definition star_rot.C:310
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Definition star_rot.h:237
virtual void equilibrium(double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff, Param *=0x0)
Computes an equilibrium configuration.
const Scalar & get_nuf() const
Returns the part of the Metric potential = logn generated by the matter terms.
Definition star_rot.h:352
const Scalar & get_bbb() const
Returns the metric factor B.
Definition star_rot.h:330
double * p_r_circ
Circumferential radius (equator).
Definition star_rot.h:248
double * p_r_circ_merid
Circumferential radius (meridian).
Definition star_rot.h:249
virtual double f_eq() const
Orbital frequency at the equator.
virtual double grv2() const
Error on the virial identity GRV2.
virtual double angu_mom() const
Angular momentum.
const Scalar & get_tnphi() const
Returns the component of the shift vector.
Definition star_rot.h:344
double * p_z_eqb
Backward redshift factor at equator.
Definition star_rot.h:253
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition star_rot.C:613
double * p_z_pole
Redshift factor at North pole.
Definition star_rot.h:254
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition star_rot.h:172
static double lambda_grv2(const Scalar &sou_m, const Scalar &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Definition star_rot.C:779
Base class for stars.
Definition star.h:175
Class intended to describe valence-2 symmetric tensors.
Definition sym_tensor.h:226
Basic array class.
Definition tbl.h:161
Tensor field of valence 1.
Definition vector.h:188
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142