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.11 2026/01/20 08:35:28 j_novak Exp $
33  * $Log: star_rot.h,v $
34  * Revision 1.11 2026/01/20 08:35:28 j_novak
35  * New method radius_Morsink() to compute the physical radius along the definition of Morsink et al. (2007) in the oblate Schwarschild approximation.
36  *
37  * Revision 1.10 2023/07/04 09:01:09 j_novak
38  * Changed the name r_pol_circ -> r_circ_merid to explicitly show the radius is a circumferential one, along a meridian.
39  *
40  * Revision 1.9 2023/07/04 07:28:34 j_novak
41  * Added method "r_pol_circ()", to compute polar circumferential radius, with a circumference along a meridian.
42  *
43  * Revision 1.8 2017/10/20 13:55:41 j_novak
44  * Calss Star_rot is now friend of class Star.
45  *
46  * Revision 1.7 2017/04/11 10:46:55 m_bejger
47  * Star_rot::surf_grav() - surface gravity values along the theta direction
48  *
49  * Revision 1.6 2015/05/19 09:30:55 j_novak
50  * New methods for computing the area of the star and its mean radius.
51  *
52  * Revision 1.5 2014/10/13 08:52:36 j_novak
53  * Lorene classes and functions now belong to the namespace Lorene.
54  *
55  * Revision 1.4 2010/02/08 10:56:30 j_novak
56  * Added a few things missing for the reading from resulting file.
57  *
58  * Revision 1.3 2010/02/02 13:35:00 e_gourgoulhon
59  * Remove the "under construction" mark.
60  *
61  * Revision 1.2 2010/01/25 18:14:05 e_gourgoulhon
62  * Added member unsurc2
63  * Added method is_relativistic()
64  * Suppressed method f_eccentric
65  *
66  * Revision 1.1 2010/01/24 16:07:45 e_gourgoulhon
67  * New class Star_rot.
68  *
69  *
70  * $Header: /cvsroot/Lorene/C++/Include/star_rot.h,v 1.11 2026/01/20 08:35:28 j_novak Exp $
71  *
72  */
73 
74 // Headers Lorene
75 #include "star.h"
76 
77 namespace Lorene {
78 class Eos ;
79 
80  //--------------------------//
81  // class Star_rot //
82  //--------------------------//
83 
101  class Star_rot : public Star {
102 
103  // Data :
104  // -----
105  protected:
110 
114  double unsurc2 ;
115 
116  double omega ;
117 
120 
123 
126 
129 
134 
137 
142 
147 
150 
153 
166 
176 
183 
202 
208 
214 
219 
224 
232 
241 
242 
243  // Derived data :
244  // ------------
245  protected:
246 
247  mutable double* p_angu_mom ;
248  mutable double* p_tsw ;
249  mutable double* p_grv2 ;
250  mutable double* p_grv3 ;
251  mutable double* p_r_circ ;
252  mutable double* p_r_circ_merid ;
253  mutable double* p_aplat ;
254  mutable double* p_area ;
255  mutable double* p_z_eqf ;
256  mutable double* p_z_eqb ;
257  mutable double* p_z_pole ;
258  mutable double* p_mom_quad ;
259  mutable double* p_r_isco ;
260  mutable double* p_f_isco ;
261  mutable double* p_espec_isco ;
264  mutable double* p_lspec_isco ;
265  mutable double* p_f_eq ;
266  mutable Tbl* p_surf_grav ;
269  mutable Tbl* p_radius_Morsink ;
270 
271 
272  // Constructors - Destructor
273  // -------------------------
274  public:
283  Star_rot(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
284 
285 
286  Star_rot(const Star_rot& ) ;
287 
295  Star_rot(Map& mp_i, const Eos& eos_i, FILE* fich) ;
296 
297  virtual ~Star_rot() ;
298 
299 
300  // Memory management
301  // -----------------
302  protected:
304  virtual void del_deriv() const ;
305 
307  virtual void set_der_0x0() const ;
308 
312  virtual void del_hydro_euler() ;
313 
314 
315  // Mutators / assignment
316  // ---------------------
317  public:
319  void operator=(const Star_rot& ) ;
320 
321  // Accessors
322  // ---------
323  public:
327  bool is_relativistic() const {return relativistic; } ;
328 
332  virtual double get_omega_c() const ;
333 
335  const Scalar& get_bbb() const {return bbb;} ;
336 
338  const Scalar& get_a_car() const {return a_car;} ;
339 
341  const Scalar& get_b_car() const {return b_car;} ;
342 
344  const Scalar& get_nphi() const {return nphi;} ;
345 
349  const Scalar& get_tnphi() const {return tnphi;} ;
350 
352  const Scalar& get_uuu() const {return uuu;} ;
353 
357  const Scalar& get_nuf() const {return nuf;} ;
358 
362  const Scalar& get_nuq() const {return nuq;} ;
363 
365  const Scalar& get_dzeta() const {return dzeta;} ;
366 
368  const Scalar& get_tggg() const {return tggg;} ;
369 
382  const Vector& get_w_shift() const {return w_shift;} ;
383 
396  const Scalar& get_khi_shift() const {return khi_shift;} ;
397 
403  const Sym_tensor& get_tkij() const {return tkij;} ;
404 
422  const Scalar& get_ak_car() const {return ak_car;} ;
423 
424  // Outputs
425  // -------
426  public:
427  virtual void sauve(FILE* ) const ;
428 
430  virtual void display_poly(ostream& ) const ;
431 
432  protected:
434  virtual ostream& operator>>(ostream& ) const ;
435 
437  virtual void partial_display(ostream& ) const ;
438 
439  // Global quantities
440  // -----------------
441  public:
442 
450  virtual const Itbl& l_surf() const ;
451 
452  virtual double mass_b() const ;
453  virtual double mass_g() const ;
454  virtual double angu_mom() const ;
455  virtual double tsw() const ;
456 
460  virtual double grv2() const ;
461 
473  virtual double grv3(ostream* ost = 0x0) const ;
474 
475  virtual double r_circ() const ;
476  virtual double r_circ_merid() const ;
477 
478  virtual double aplat() const ;
479  virtual double area() const ;
480  virtual double mean_radius() const ;
482  virtual double z_eqf() const ;
483  virtual double z_eqb() const ;
484  virtual double z_pole() const ;
485 
495  virtual double mom_quad() const ;
496 
500  virtual const Tbl& surf_grav() const ;
501 
505  virtual const Tbl& radius_Morsink() const ;
506 
513  virtual double r_isco(ostream* ost = 0x0) const ;
514 
516  virtual double f_isco() const ;
517 
519  virtual double espec_isco() const ;
520 
522  virtual double lspec_isco() const ;
523 
525  virtual double f_eq() const ;
526 
527 
528  // Computational routines
529  // ----------------------
530  public:
541  virtual void hydro_euler() ;
542 
553  void update_metric() ;
554 
563  void fait_shift() ;
564 
568  void fait_nphi() ;
569 
573  void extrinsic_curvature() ;
574 
604  static double lambda_grv2(const Scalar& sou_m, const Scalar& sou_q) ;
605 
685  virtual void equilibrium(double ent_c, double omega0, double fact_omega,
686  int nzadapt, const Tbl& ent_limit,
687  const Itbl& icontrol, const Tbl& control,
688  double mbar_wanted, double aexp_mass,
689  Tbl& diff, Param* = 0x0) ;
690 
691 
692  friend class Star ;
693  };
694 
695 }
696 #endif
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.
virtual double mass_g() const
Gravitational mass.
const Scalar & get_a_car() const
Returns the square of the metric factor A.
Definition: star_rot.h:338
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: star_rot.h:240
Scalar dzeta
Metric potential .
Definition: star_rot.h:149
virtual double mean_radius() const
Mean star radius from the area .
Scalar ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta ...
Definition: star_rot.h:218
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:119
virtual double lspec_isco() const
Angular momentum of a particle on the ISCO.
static double lambda_grv2(const Scalar &sou_m, const Scalar &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition: star_rot.h:264
double * p_r_isco
Circumferential radius of the ISCO.
Definition: star_rot.h:259
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: star_rot.h:141
virtual double r_circ_merid() const
Crcumferential meridional radius.
double * p_f_eq
Orbital frequency at the equator.
Definition: star_rot.h:265
double * p_r_circ
Circumferential radius (equator)
Definition: star_rot.h:251
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition: star_rot.h:175
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: star_rot.h:109
Scalar tnphi
Component of the shift vector.
Definition: star_rot.h:133
double * p_mom_quad
Quadrupole moment.
Definition: star_rot.h:258
virtual double f_eq() const
Orbital frequency at the equator.
Scalar bbb
Metric factor B.
Definition: star_rot.h:122
Lorene prototypes.
Definition: app_hor.h:67
Equation of state base class.
Definition: eos.h:209
virtual double z_pole() const
Redshift factor at North pole.
double * p_f_isco
Orbital frequency of the ISCO.
Definition: star_rot.h:260
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
virtual double grv2() const
Error on the virial identity GRV2.
Base class for coordinate mappings.
Definition: map.h:697
void extrinsic_curvature()
Computes tkij and ak_car from shift , nnn and b_car .
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:671
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star_rot.C:313
Basic integer array class.
Definition: itbl.h:122
void update_metric()
Computes metric coefficients from known potentials.
Scalar nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: star_rot.h:146
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot.h:116
Base class for stars.
Definition: star.h:175
virtual const Tbl & surf_grav() const
Surface gravity (table along the theta direction)
Tensor field of valence 1.
Definition: vector.h:188
Tbl * p_radius_Morsink
Physical radius close to the definition of Morsink et al. (2007)
Definition: star_rot.h:269
bool is_relativistic() const
Returns true for a relativistic star, false for a Newtonian one.
Definition: star_rot.h:327
virtual double area() const
Integrated surface area in .
Class for isolated rotating stars.
Definition: star_rot.h:101
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...
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Definition: star_rot.h:223
const Vector & get_w_shift() const
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog...
Definition: star_rot.h:382
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: star_rot.C:341
void operator=(const Star_rot &)
Assignment to another Star_rot.
Definition: star_rot.C:381
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star_rot.C:450
virtual ~Star_rot()
Destructor.
Definition: star_rot.C:303
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
Definition: star_rot_isco.C:74
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition: star_rot.h:114
double * p_tsw
Ratio T/W.
Definition: star_rot.h:248
Scalar nphi
Metric coefficient .
Definition: star_rot.h:128
virtual double espec_isco() const
Energy of a particle on the ISCO.
const Scalar & get_ak_car() const
Returns the scalar .
Definition: star_rot.h:422
double * p_z_eqf
Forward redshift factor at equator.
Definition: star_rot.h:255
const Scalar & get_khi_shift() const
Returns the scalar used in the decomposition of shift following Shibata&#39;s prescription [Prog...
Definition: star_rot.h:396
double * p_grv2
Error on the virial identity GRV2.
Definition: star_rot.h:249
double * p_area
Integrated surface area.
Definition: star_rot.h:254
double * p_r_circ_merid
Circumferential radius (meridian)
Definition: star_rot.h:252
virtual double mass_b() const
Baryon mass.
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition: star_rot.h:262
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: star_rot.C:367
virtual double f_isco() const
Orbital frequency at the innermost stable circular orbit (ISCO).
virtual const Tbl & radius_Morsink() const
Physical radius following Morsink et al.
virtual double mom_quad() const
Quadrupole moment.
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:213
virtual double angu_mom() const
Angular momentum.
double * p_z_pole
Redshift factor at North pole.
Definition: star_rot.h:257
Parameter storage.
Definition: param.h:125
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: star_rot.C:680
Star_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition: star_rot.C:100
const Scalar & get_bbb() const
Returns the metric factor B.
Definition: star_rot.h:335
double * p_grv3
Error on the virial identity GRV3.
Definition: star_rot.h:250
const Scalar & get_dzeta() const
Returns the Metric potential .
Definition: star_rot.h:365
Sym_tensor tkij
Tensor related to the extrinsic curvature tensor by .
Definition: star_rot.h:182
const Sym_tensor & get_tkij() const
Returns the tensor related to the extrinsic curvature tensor by .
Definition: star_rot.h:403
virtual double r_circ() const
Circumferential equatorial radius.
double * p_angu_mom
Angular momentum.
Definition: star_rot.h:247
Scalar b_car
Square of the metric factor B.
Definition: star_rot.h:125
const Scalar & get_tnphi() const
Returns the component of the shift vector.
Definition: star_rot.h:349
const Scalar & get_nuq() const
Returns the Part of the Metric potential = logn generated by the quadratic terms.
Definition: star_rot.h:362
const Scalar & get_nuf() const
Returns the part of the Metric potential = logn generated by the matter terms.
Definition: star_rot.h:357
double * p_aplat
Flatening r_pole/r_eq.
Definition: star_rot.h:253
virtual double z_eqb() const
Backward redshift factor at equator.
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition: star_rot.C:618
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:231
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:136
const Scalar & get_nphi() const
Returns the metric coefficient .
Definition: star_rot.h:344
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:207
virtual void sauve(FILE *) const
Save in a file.
Definition: star_rot.C:423
Basic array class.
Definition: tbl.h:164
const Scalar & get_uuu() const
Returns the norm of u_euler.
Definition: star_rot.h:352
Scalar tggg
Metric potential .
Definition: star_rot.h:152
const Scalar & get_b_car() const
Returns the square of the metric factor B.
Definition: star_rot.h:341
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Definition: star_rot.C:784
Scalar ak_car
Scalar .
Definition: star_rot.h:201
virtual double aplat() const
Flatening r_pole/r_eq.
Tbl * p_surf_grav
Surface gravity (along the theta direction)
Definition: star_rot.h:267
Vector w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: star_rot.h:165
const Scalar & get_tggg() const
Returns the Metric potential .
Definition: star_rot.h:368
virtual double z_eqf() const
Forward redshift factor at equator.
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata&#39;s prescription [Prog.
Definition: star_rot.C:735
virtual double tsw() const
Ratio T/W.
double * p_z_eqb
Backward redshift factor at equator.
Definition: star_rot.h:256