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.12 2026/02/05 14:23:25 j_novak Exp $
33  * $Log: star_rot.h,v $
34  * Revision 1.12 2026/02/05 14:23:25 j_novak
35  * New method surf_grav_euler().
36  *
37  * Revision 1.11 2026/01/20 08:35:28 j_novak
38  * New method radius_Morsink() to compute the physical radius along the definition of Morsink et al. (2007) in the oblate Schwarschild approximation.
39  *
40  * Revision 1.10 2023/07/04 09:01:09 j_novak
41  * Changed the name r_pol_circ -> r_circ_merid to explicitly show the radius is a circumferential one, along a meridian.
42  *
43  * Revision 1.9 2023/07/04 07:28:34 j_novak
44  * Added method "r_pol_circ()", to compute polar circumferential radius, with a circumference along a meridian.
45  *
46  * Revision 1.8 2017/10/20 13:55:41 j_novak
47  * Calss Star_rot is now friend of class Star.
48  *
49  * Revision 1.7 2017/04/11 10:46:55 m_bejger
50  * Star_rot::surf_grav() - surface gravity values along the theta direction
51  *
52  * Revision 1.6 2015/05/19 09:30:55 j_novak
53  * New methods for computing the area of the star and its mean radius.
54  *
55  * Revision 1.5 2014/10/13 08:52:36 j_novak
56  * Lorene classes and functions now belong to the namespace Lorene.
57  *
58  * Revision 1.4 2010/02/08 10:56:30 j_novak
59  * Added a few things missing for the reading from resulting file.
60  *
61  * Revision 1.3 2010/02/02 13:35:00 e_gourgoulhon
62  * Remove the "under construction" mark.
63  *
64  * Revision 1.2 2010/01/25 18:14:05 e_gourgoulhon
65  * Added member unsurc2
66  * Added method is_relativistic()
67  * Suppressed method f_eccentric
68  *
69  * Revision 1.1 2010/01/24 16:07:45 e_gourgoulhon
70  * New class Star_rot.
71  *
72  *
73  * $Header: /cvsroot/Lorene/C++/Include/star_rot.h,v 1.12 2026/02/05 14:23:25 j_novak Exp $
74  *
75  */
76 
77 // Headers Lorene
78 #include "star.h"
79 
80 namespace Lorene {
81 class Eos ;
82 
83  //--------------------------//
84  // class Star_rot //
85  //--------------------------//
86 
104  class Star_rot : public Star {
105 
106  // Data :
107  // -----
108  protected:
113 
117  double unsurc2 ;
118 
119  double omega ;
120 
123 
126 
129 
132 
137 
140 
145 
150 
153 
156 
169 
179 
186 
205 
211 
217 
222 
227 
235 
244 
245 
246  // Derived data :
247  // ------------
248  protected:
249 
250  mutable double* p_angu_mom ;
251  mutable double* p_tsw ;
252  mutable double* p_grv2 ;
253  mutable double* p_grv3 ;
254  mutable double* p_r_circ ;
255  mutable double* p_r_circ_merid ;
256  mutable double* p_aplat ;
257  mutable double* p_area ;
258  mutable double* p_z_eqf ;
259  mutable double* p_z_eqb ;
260  mutable double* p_z_pole ;
261  mutable double* p_mom_quad ;
262  mutable double* p_r_isco ;
263  mutable double* p_f_isco ;
264  mutable double* p_espec_isco ;
267  mutable double* p_lspec_isco ;
268  mutable double* p_f_eq ;
269  mutable Tbl* p_surf_grav ;
274  mutable Tbl* p_radius_Morsink ;
275 
276 
277  // Constructors - Destructor
278  // -------------------------
279  public:
288  Star_rot(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
289 
290 
291  Star_rot(const Star_rot& ) ;
292 
300  Star_rot(Map& mp_i, const Eos& eos_i, FILE* fich) ;
301 
302  virtual ~Star_rot() ;
303 
304 
305  // Memory management
306  // -----------------
307  protected:
309  virtual void del_deriv() const ;
310 
312  virtual void set_der_0x0() const ;
313 
317  virtual void del_hydro_euler() ;
318 
319 
320  // Mutators / assignment
321  // ---------------------
322  public:
324  void operator=(const Star_rot& ) ;
325 
326  // Accessors
327  // ---------
328  public:
332  bool is_relativistic() const {return relativistic; } ;
333 
337  virtual double get_omega_c() const ;
338 
340  const Scalar& get_bbb() const {return bbb;} ;
341 
343  const Scalar& get_a_car() const {return a_car;} ;
344 
346  const Scalar& get_b_car() const {return b_car;} ;
347 
349  const Scalar& get_nphi() const {return nphi;} ;
350 
354  const Scalar& get_tnphi() const {return tnphi;} ;
355 
357  const Scalar& get_uuu() const {return uuu;} ;
358 
362  const Scalar& get_nuf() const {return nuf;} ;
363 
367  const Scalar& get_nuq() const {return nuq;} ;
368 
370  const Scalar& get_dzeta() const {return dzeta;} ;
371 
373  const Scalar& get_tggg() const {return tggg;} ;
374 
387  const Vector& get_w_shift() const {return w_shift;} ;
388 
401  const Scalar& get_khi_shift() const {return khi_shift;} ;
402 
408  const Sym_tensor& get_tkij() const {return tkij;} ;
409 
427  const Scalar& get_ak_car() const {return ak_car;} ;
428 
429  // Outputs
430  // -------
431  public:
432  virtual void sauve(FILE* ) const ;
433 
435  virtual void display_poly(ostream& ) const ;
436 
437  protected:
439  virtual ostream& operator>>(ostream& ) const ;
440 
442  virtual void partial_display(ostream& ) const ;
443 
444  // Global quantities
445  // -----------------
446  public:
447 
455  virtual const Itbl& l_surf() const ;
456 
457  virtual double mass_b() const ;
458  virtual double mass_g() const ;
459  virtual double angu_mom() const ;
460  virtual double tsw() const ;
461 
465  virtual double grv2() const ;
466 
478  virtual double grv3(ostream* ost = 0x0) const ;
479 
480  virtual double r_circ() const ;
481  virtual double r_circ_merid() const ;
482 
483  virtual double aplat() const ;
484  virtual double area() const ;
485  virtual double mean_radius() const ;
487  virtual double z_eqf() const ;
488  virtual double z_eqb() const ;
489  virtual double z_pole() const ;
490 
500  virtual double mom_quad() const ;
501 
505  virtual const Tbl& surf_grav() const ;
506 
510  virtual const Tbl& surf_grav_euler() const ;
511 
515  virtual const Tbl& radius_Morsink() const ;
516 
523  virtual double r_isco(ostream* ost = 0x0) const ;
524 
526  virtual double f_isco() const ;
527 
529  virtual double espec_isco() const ;
530 
532  virtual double lspec_isco() const ;
533 
535  virtual double f_eq() const ;
536 
537 
538  // Computational routines
539  // ----------------------
540  public:
551  virtual void hydro_euler() ;
552 
563  void update_metric() ;
564 
573  void fait_shift() ;
574 
578  void fait_nphi() ;
579 
583  void extrinsic_curvature() ;
584 
614  static double lambda_grv2(const Scalar& sou_m, const Scalar& sou_q) ;
615 
695  virtual void equilibrium(double ent_c, double omega0, double fact_omega,
696  int nzadapt, const Tbl& ent_limit,
697  const Itbl& icontrol, const Tbl& control,
698  double mbar_wanted, double aexp_mass,
699  Tbl& diff, Param* = 0x0) ;
700 
701 
702  friend class Star ;
703  };
704 
705 }
706 #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:343
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: star_rot.h:243
Scalar dzeta
Metric potential .
Definition: star_rot.h:152
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:221
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:122
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:267
double * p_r_isco
Circumferential radius of the ISCO.
Definition: star_rot.h:262
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: star_rot.h:144
virtual double r_circ_merid() const
Crcumferential meridional radius.
double * p_f_eq
Orbital frequency at the equator.
Definition: star_rot.h:268
double * p_r_circ
Circumferential radius (equator)
Definition: star_rot.h:254
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition: star_rot.h:178
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: star_rot.h:112
Scalar tnphi
Component of the shift vector.
Definition: star_rot.h:136
double * p_mom_quad
Quadrupole moment.
Definition: star_rot.h:261
virtual double f_eq() const
Orbital frequency at the equator.
Scalar bbb
Metric factor B.
Definition: star_rot.h:125
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:263
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 const Tbl & surf_grav_euler() const
Surface gravity for the Eulerian observer (table along the theta direction)
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:676
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star_rot.C:316
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:149
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot.h:119
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:274
bool is_relativistic() const
Returns true for a relativistic star, false for a Newtonian one.
Definition: star_rot.h:332
virtual double area() const
Integrated surface area in .
Class for isolated rotating stars.
Definition: star_rot.h:104
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:226
const Vector & get_w_shift() const
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog...
Definition: star_rot.h:387
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: star_rot.C:345
void operator=(const Star_rot &)
Assignment to another Star_rot.
Definition: star_rot.C:386
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star_rot.C:455
virtual ~Star_rot()
Destructor.
Definition: star_rot.C:306
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:117
double * p_tsw
Ratio T/W.
Definition: star_rot.h:251
Scalar nphi
Metric coefficient .
Definition: star_rot.h:131
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:427
double * p_z_eqf
Forward redshift factor at equator.
Definition: star_rot.h:258
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:401
double * p_grv2
Error on the virial identity GRV2.
Definition: star_rot.h:252
double * p_area
Integrated surface area.
Definition: star_rot.h:257
double * p_r_circ_merid
Circumferential radius (meridian)
Definition: star_rot.h:255
virtual double mass_b() const
Baryon mass.
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition: star_rot.h:265
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: star_rot.C:372
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:216
virtual double angu_mom() const
Angular momentum.
double * p_z_pole
Redshift factor at North pole.
Definition: star_rot.h:260
Parameter storage.
Definition: param.h:125
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: star_rot.C:685
Star_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition: star_rot.C:103
const Scalar & get_bbb() const
Returns the metric factor B.
Definition: star_rot.h:340
double * p_grv3
Error on the virial identity GRV3.
Definition: star_rot.h:253
const Scalar & get_dzeta() const
Returns the Metric potential .
Definition: star_rot.h:370
Tbl * p_surf_grav_euler
Surface gravity for the Eulerian observer (along the theta direction)
Definition: star_rot.h:272
Sym_tensor tkij
Tensor related to the extrinsic curvature tensor by .
Definition: star_rot.h:185
const Sym_tensor & get_tkij() const
Returns the tensor related to the extrinsic curvature tensor by .
Definition: star_rot.h:408
virtual double r_circ() const
Circumferential equatorial radius.
double * p_angu_mom
Angular momentum.
Definition: star_rot.h:250
Scalar b_car
Square of the metric factor B.
Definition: star_rot.h:128
const Scalar & get_tnphi() const
Returns the component of the shift vector.
Definition: star_rot.h:354
const Scalar & get_nuq() const
Returns the Part of the Metric potential = logn generated by the quadratic terms.
Definition: star_rot.h:367
const Scalar & get_nuf() const
Returns the part of the Metric potential = logn generated by the matter terms.
Definition: star_rot.h:362
double * p_aplat
Flatening r_pole/r_eq.
Definition: star_rot.h:256
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:623
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:234
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:139
const Scalar & get_nphi() const
Returns the metric coefficient .
Definition: star_rot.h:349
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:210
virtual void sauve(FILE *) const
Save in a file.
Definition: star_rot.C:428
Basic array class.
Definition: tbl.h:164
const Scalar & get_uuu() const
Returns the norm of u_euler.
Definition: star_rot.h:357
Scalar tggg
Metric potential .
Definition: star_rot.h:155
const Scalar & get_b_car() const
Returns the square of the metric factor B.
Definition: star_rot.h:346
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Definition: star_rot.C:789
Scalar ak_car
Scalar .
Definition: star_rot.h:204
virtual double aplat() const
Flatening r_pole/r_eq.
Tbl * p_surf_grav
Surface gravity (along the theta direction)
Definition: star_rot.h:270
Vector w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: star_rot.h:168
const Scalar & get_tggg() const
Returns the Metric potential .
Definition: star_rot.h:373
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:740
virtual double tsw() const
Ratio T/W.
double * p_z_eqb
Backward redshift factor at equator.
Definition: star_rot.h:259