LORENE
star_rot_global.C
1 /*
2  * Methods for computing global quantities within the class Star_rot
3  *
4  * (see file star_rot.h for documentation)
5  */
6 
7 /*
8  * Copyright (c) 2010 Eric Gourgoulhon
9  *
10  * This file is part of LORENE.
11  *
12  * LORENE is free software; you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation; either version 2 of the License, or
15  * (at your option) any later version.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 
30 
31 /*
32  * $Id: star_rot_global.C,v 1.13 2026/01/20 08:35:28 j_novak Exp $
33  * $Log: star_rot_global.C,v $
34  * Revision 1.13 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.12 2026/01/19 13:04:23 j_novak
38  * Corrected the location of the surface in surf_gravity()
39  *
40  * Revision 1.11 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.10 2023/07/03 14:32:31 j_novak
44  * Added method "r_pol_circ()", to compute polar circumferential radius, with a circumference along a meridian.
45  *
46  * Revision 1.9 2021/10/06 09:41:28 j_novak
47  * Corrected quadrupole computation from Friedman & Stergioulas book
48  *
49  * Revision 1.8 2017/04/13 09:59:34 m_bejger
50  * Bug fix (contraction using physical gamma metric in surf_grav - vector a)
51  *
52  * Revision 1.7 2017/04/11 10:46:55 m_bejger
53  * Star_rot::surf_grav() - surface gravity values along the theta direction
54  *
55  * Revision 1.6 2016/12/05 16:18:15 j_novak
56  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
57  *
58  * Revision 1.5 2015/05/19 09:30:56 j_novak
59  * New methods for computing the area of the star and its mean radius.
60  *
61  * Revision 1.4 2014/10/13 08:53:39 j_novak
62  * Lorene classes and functions now belong to the namespace Lorene.
63  *
64  * Revision 1.3 2014/10/06 15:13:17 j_novak
65  * Modified #include directives to use c++ syntax.
66  *
67  * Revision 1.2 2010/01/25 22:33:35 e_gourgoulhon
68  * Debugging...
69  *
70  * Revision 1.1 2010/01/25 18:15:52 e_gourgoulhon
71  * First version.
72  *
73  *
74  *
75  * $Header: /cvsroot/Lorene/C++/Source/Star/star_rot_global.C,v 1.13 2026/01/20 08:35:28 j_novak Exp $
76  *
77  */
78 
79 // Headers C
80 #include <cstdlib>
81 #include <cmath>
82 
83 // Headers Lorene
84 #include "star_rot.h"
85 #include "unites.h"
86 
87  //--------------------------//
88  // Stellar surface //
89  //--------------------------//
90 
91 namespace Lorene {
92  const Itbl& Star_rot::l_surf() const {
93 
94  if (p_l_surf == 0x0) { // a new computation is required
95 
96  assert(p_xi_surf == 0x0) ; // consistency check
97 
98  int np = mp.get_mg()->get_np(0) ;
99  int nt = mp.get_mg()->get_nt(0) ;
100 
101  p_l_surf = new Itbl(np, nt) ;
102  p_xi_surf = new Tbl(np, nt) ;
103 
104  double ent0 = 0 ; // definition of the surface
105  double precis = 1.e-15 ;
106  int nitermax = 100 ;
107  int niter ;
108 
109  (ent.get_spectral_va()).equipot(ent0, nzet, precis, nitermax, niter, *p_l_surf,
110  *p_xi_surf) ;
111 
112  }
113 
114  return *p_l_surf ;
115 
116  }
117 
118  //------------------------------//
119  // Baryon mass //
120  //------------------------------//
121 
122  double Star_rot::mass_b() const {
123 
124  if (p_mass_b == 0x0) { // a new computation is required
125 
126  if (relativistic) {
127 
128  Scalar dens = a_car * bbb * gam_euler * nbar ;
129 
130  p_mass_b = new double( dens.integrale() ) ;
131  }
132  else{ // Newtonian case
133  assert(nbar.get_etat() == ETATQCQ) ;
134 
135  p_mass_b = new double( nbar.integrale() ) ;
136 
137  }
138 
139  }
140  return *p_mass_b ;
141  }
142 
143 
144  //----------------------------//
145  // Gravitational mass //
146  //----------------------------//
147 
148  double Star_rot::mass_g() const {
149 
150  if (p_mass_g == 0x0) { // a new computation is required
151 
152  if (relativistic) {
153 
154  Scalar source = nn * (ener_euler + s_euler)
155  + 2 * bbb * (ener_euler + press)
156  * tnphi * uuu ;
157  source = a_car * bbb * source ;
158  source.std_spectral_base() ;
159 
160  p_mass_g = new double( source.integrale() ) ;
161  }
162  else{ // Newtonian case
163  p_mass_g = new double( mass_b() ) ; // in the Newtonian case
164  // M_g = M_b
165  }
166  }
167  return *p_mass_g ;
168 }
169 
170  //----------------------------//
171  // Angular momentum //
172  //----------------------------//
173 
174  double Star_rot::angu_mom() const {
175 
176  if (p_angu_mom == 0x0) { // a new computation is required
177 
178  Scalar dens = uuu ;
179 
180  dens.mult_r() ; // Multiplication by
181  dens.set_spectral_va() = (dens.get_spectral_va()).mult_st() ; // r sin(theta)
182 
183  if (relativistic) {
184  dens = a_car * b_car * (ener_euler + press) * dens ;
185  }
186  else { // Newtonian case
187  dens = nbar * dens ;
188  }
189 
190  p_angu_mom = new double( dens.integrale() ) ;
191  }
192  return *p_angu_mom ;
193  }
194 
195 
196  //----------------------------//
197  // T/W //
198  //----------------------------//
199 
200  double Star_rot::tsw() const {
201 
202  if (p_tsw == 0x0) { // a new computation is required
203 
204  double tcin = 0.5 * omega * angu_mom() ;
205 
206  if (relativistic) {
207 
208  Scalar dens = a_car * bbb * gam_euler * ener ;
209  double mass_p = dens.integrale() ;
210 
211  p_tsw = new double( tcin / ( mass_p + tcin - mass_g() ) ) ;
212 
213  }
214  else { // Newtonian case
215  Scalar dens = 0.5 * nbar * logn ;
216  double wgrav = dens.integrale() ;
217  p_tsw = new double( tcin / fabs(wgrav) ) ;
218  }
219  }
220  return *p_tsw ;
221 }
222 
223 
224  //----------------------------//
225  // GRV2 //
226  //----------------------------//
227 
228 double Star_rot::grv2() const {
229 
230  using namespace Unites ;
231  if (p_grv2 == 0x0) { // a new computation is required
232 
233  Scalar sou_m(mp) ;
234  if (relativistic) {
235  sou_m = 2 * qpig * a_car * (press + (ener_euler+press)
236  * uuu*uuu ) ;
237  }
238  else {
239  sou_m = 2 * qpig * (press + nbar * uuu*uuu ) ;
240  }
241 
242  Vector dlogn = logn.derive_cov( mp.flat_met_spher() ) ;
243  Scalar sou_q = 1.5 * ak_car
244  - dlogn(1)*dlogn(1) - dlogn(2)*dlogn(2) - dlogn(3)*dlogn(3) ;
245 
246  p_grv2 = new double( double(1) - lambda_grv2(sou_m, sou_q) ) ;
247 
248  }
249 
250  return *p_grv2 ;
251 
252 }
253 
254 
255  //----------------------------//
256  // GRV3 //
257  //----------------------------//
258 
259 double Star_rot::grv3(ostream* ost) const {
260 
261  using namespace Unites ;
262 
263  if (p_grv3 == 0x0) { // a new computation is required
264 
265  Scalar source(mp) ;
266 
267  // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
268  // ------------------ Class. Quantum Grav. 11, 443 (1994)]
269 
270  Vector dlogn = logn.derive_cov( mp.flat_met_spher() ) ;
271 
272  if (relativistic) {
273  Scalar alpha = dzeta - logn ;
274  Scalar bet = log( bbb ) ;
275  bet.std_spectral_base() ;
276 
277  Vector dalpha = alpha.derive_cov( mp.flat_met_spher() ) ;
278  Vector dbet = bet.derive_cov( mp.flat_met_spher() ) ;
279 
280  source = 0.75 * ak_car
281  - dlogn(1)*dlogn(1) - dlogn(2)*dlogn(2) - dlogn(3)*dlogn(3)
282  + 0.5 * ( dalpha(1)*dbet(1) + dalpha(2)*dbet(2) + dalpha(3)*dbet(3) ) ;
283 
284  Scalar aa = alpha - 0.5 * bet ;
285  Scalar daadt = aa.srdsdt() ; // 1/r d/dth
286 
287  // What follows is valid only for a mapping of class Map_radial :
288  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
289  if (mpr == 0x0) {
290  cout << "Star_rot::grv3: the mapping does not belong"
291  << " to the class Map_radial !" << endl ;
292  abort() ;
293  }
294 
295  // Computation of 1/tan(theta) * 1/r daa/dtheta
296  if (daadt.get_etat() == ETATQCQ) {
297  Valeur& vdaadt = daadt.set_spectral_va() ;
298  vdaadt = vdaadt.ssint() ; // division by sin(theta)
299  vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
300  }
301 
302  Scalar temp = aa.dsdr() + daadt ;
303  temp = ( bbb - a_car/bbb ) * temp ;
304  temp.std_spectral_base() ;
305 
306  // Division by r
307  Valeur& vtemp = temp.set_spectral_va() ;
308  vtemp = vtemp.sx() ; // division by xi in the nucleus
309  // Id in the shells
310  // division by xi-1 in the ZEC
311  vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
312  // by 1/r in the shells
313  // by r(xi-1) in the ZEC
314 
315  // In the ZEC, a multiplication by r has been performed instead
316  // of the division:
317  temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;
318 
319  source = bbb * source + 0.5 * temp ;
320 
321  }
322  else{
323  source = - 0.5 * ( dlogn(1)*dlogn(1) + dlogn(2)*dlogn(2) + dlogn(3)*dlogn(3) ) ;
324  }
325 
326  source.std_spectral_base() ;
327 
328  double int_grav = source.integrale() ;
329 
330  // Matter term
331  // -----------
332 
333  if (relativistic) {
334  source = qpig * a_car * bbb * s_euler ;
335  }
336  else{
337  source = qpig * ( 3 * press + nbar * uuu * uuu ) ;
338  }
339 
340  source.std_spectral_base() ;
341 
342  double int_mat = source.integrale() ;
343 
344  // Virial error
345  // ------------
346  if (ost != 0x0) {
347  *ost << "Star_rot::grv3 : gravitational term : " << int_grav
348  << endl ;
349  *ost << "Star_rot::grv3 : matter term : " << int_mat
350  << endl ;
351  }
352 
353  p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;
354 
355  }
356 
357  return *p_grv3 ;
358 
359 }
360 
361 
362  //----------------------------//
363  // R_circ //
364  //----------------------------//
365 
366 double Star_rot::r_circ() const {
367 
368  if (p_r_circ == 0x0) { // a new computation is required
369 
370  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
371  const Mg3d* mg = mp.get_mg() ;
372  assert(mg->get_type_t() == SYM) ;
373  int l_b = nzet - 1 ;
374  int i_b = mg->get_nr(l_b) - 1 ;
375  int j_b = mg->get_nt(l_b) - 1 ;
376  int k_b = 0 ;
377 
378  p_r_circ = new double( bbb.val_grid_point(l_b, k_b, j_b, i_b) * ray_eq() ) ;
379 
380  }
381 
382  return *p_r_circ ;
383 
384 }
385 
386 double Star_rot::r_circ_merid() const {
387 
388  if (p_r_circ_merid == 0x0) { // a new computation is required
389  const Mg3d& mg = *(mp.get_mg()) ;
390  int np = mg.get_np(0) ;
391  int nt = mg.get_nt(0) ;
392  assert(np == 1) ; //Only valid for axisymmetric configurations
393 
394  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
395  assert(mp_rad != 0x0) ;
396 
397  Valeur va_r(mg.get_angu()) ;
398  va_r.annule_hard() ;
399  Itbl lsurf = l_surf() ;
400  Tbl xisurf = xi_surf() ;
401  for (int k=0; k<np; k++) {
402  for (int j=0; j<nt; j++) {
403  int l_star = lsurf(k, j) ;
404  double xi_star = xisurf(k, j) ;
405  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star, xi_star, j, k) ;
406  }
407  }
408  va_r.std_base_scal() ;
409 
410  Valeur integ(mg.get_angu()) ;
411  Valeur dr = va_r.dsdt() ;
412  integ = sqrt(va_r*va_r + dr*dr) ;
414  int k = 0 ;
415  for (int j=0; j<nt; j++) {
416  integ.set(0, k, j, 0)
417  *= sqrt(a2.val_point_jk(lsurf(k, j), xisurf(k, j), j, k)) ;
418  }
419  integ.std_base_scal() ;
420  integ.coef() ;
421  double surftot = (*integ.c_cf)(0, k, 0, 0) ; // Only j=0 has a nonvanishing
422  // integral
423  p_r_circ_merid = new double(surftot ) ;
424 
425  }
426 
427  return *p_r_circ_merid ;
428 
429 }
430 
431  //----------------------------//
432  // Surface area //
433  //----------------------------//
434 
435  double Star_rot::area() const {
436 
437  if (p_area == 0x0) { // a new computation is required
438  const Mg3d& mg = *(mp.get_mg()) ;
439  int np = mg.get_np(0) ;
440  int nt = mg.get_nt(0) ;
441  assert(np == 1) ; //Only valid for axisymmetric configurations
442 
443  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
444  assert(mp_rad != 0x0) ;
445 
446  Valeur va_r(mg.get_angu()) ;
447  va_r.annule_hard() ;
448  Itbl lsurf = l_surf() ;
449  Tbl xisurf = xi_surf() ;
450  for (int k=0; k<np; k++) {
451  for (int j=0; j<nt; j++) {
452  int l_star = lsurf(k, j) ;
453  double xi_star = xisurf(k, j) ;
454  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star, xi_star, j, k) ;
455  }
456  }
457  va_r.std_base_scal() ;
458 
459  Valeur integ(mg.get_angu()) ;
460  Valeur dr = va_r.dsdt() ;
461  integ = sqrt(va_r*va_r + dr*dr) ;
464  for (int k=0; k<np; k++) {
465  for (int j=0; j<nt; j++) {
466  integ.set(0, k, j, 0) *= sqrt(a2.val_point_jk(lsurf(k, j), xisurf(k, j), j, k))
467  * b.val_point_jk(lsurf(k, j), xisurf(k, j), j, k) * va_r(0, k, j, 0) ;
468  }
469  }
470  integ.std_base_scal() ;
471  Valeur integ2 = integ.mult_st() ;
472  double surftot = 0. ;
473  for (int j=0; j<nt; j++) {
474  surftot += (*integ2.c_cf)(0, 0, j, 0) / double(2*j+1) ;
475  }
476 
477  p_area = new double( 4*M_PI*surftot) ;
478 
479  }
480 
481  return *p_area ;
482 
483 }
484 
485  double Star_rot::mean_radius() const {
486 
487  return sqrt(area()/(4*M_PI)) ;
488 
489  }
490 
491  //----------------------------//
492  // Flattening //
493  //----------------------------//
494 
495 double Star_rot::aplat() const {
496 
497  if (p_aplat == 0x0) { // a new computation is required
498 
499  p_aplat = new double( ray_pole() / ray_eq() ) ;
500 
501  }
502 
503  return *p_aplat ;
504 
505 }
506 
507 
508  //----------------------------//
509  // Z_eq_f //
510  //----------------------------//
511 
512 double Star_rot::z_eqf() const {
513 
514  if (p_z_eqf == 0x0) { // a new computation is required
515 
516  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
517  const Mg3d* mg = mp.get_mg() ;
518  assert(mg->get_type_t() == SYM) ;
519  int l_b = nzet - 1 ;
520  int i_b = mg->get_nr(l_b) - 1 ;
521  int j_b = mg->get_nt(l_b) - 1 ;
522  int k_b = 0 ;
523 
524  double u_eq = uuu.val_grid_point(l_b, k_b, j_b, i_b) ;
525  double n_eq = nn.val_grid_point(l_b, k_b, j_b, i_b) ;
526  double nphi_eq = nphi.val_grid_point(l_b, k_b, j_b, i_b) ;
527 
528  p_z_eqf = new double( sqrt((1.-u_eq)/(1.+u_eq))
529  / (n_eq + nphi_eq * r_circ() )
530  - 1. ) ;
531  }
532 
533  return *p_z_eqf ;
534 
535 }
536  //----------------------------//
537  // Z_eq_b //
538  //----------------------------//
539 
540 double Star_rot::z_eqb() const {
541 
542  if (p_z_eqb == 0x0) { // a new computation is required
543 
544  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
545  const Mg3d* mg = mp.get_mg() ;
546  assert(mg->get_type_t() == SYM) ;
547  int l_b = nzet - 1 ;
548  int i_b = mg->get_nr(l_b) - 1 ;
549  int j_b = mg->get_nt(l_b) - 1 ;
550  int k_b = 0 ;
551 
552  double u_eq = uuu.val_grid_point(l_b, k_b, j_b, i_b) ;
553  double n_eq = nn.val_grid_point(l_b, k_b, j_b, i_b) ;
554  double nphi_eq = nphi.val_grid_point(l_b, k_b, j_b, i_b) ;
555 
556  p_z_eqb = new double( sqrt((1.+u_eq)/(1.-u_eq))
557  / (n_eq - nphi_eq * r_circ() )
558  - 1. ) ;
559 
560  }
561 
562  return *p_z_eqb ;
563 
564 }
565 
566 
567  //----------------------------//
568  // Z_pole //
569  //----------------------------//
570 
571 double Star_rot::z_pole() const {
572 
573  if (p_z_pole == 0x0) { // a new computation is required
574 
575  double n_pole = nn.val_point(ray_pole(), 0., 0.) ;
576 
577  p_z_pole = new double( 1. / n_pole - 1. ) ;
578 
579  }
580 
581  return *p_z_pole ;
582 
583 }
584 
585 
586  //----------------------------//
587  // Quadrupole moment //
588  // Uses Eq.(11) of Pappas & //
589  // Apostolatos, PRL 108 2012 //
590  //----------------------------//
591 
592 double Star_rot::mom_quad() const {
593 
594  using namespace Unites ;
595  if (p_mom_quad == 0x0) { // a new computation is required
596 
597  // Source for of the Poisson equation for nu
598  // -----------------------------------------
599 
600  Scalar source(mp) ;
601 
602  if (relativistic) {
603  Scalar bet = log(bbb) ;
604  bet.std_spectral_base() ;
605 
606  Vector dlogn = logn.derive_cov( mp.flat_met_spher() ) ;
607  Vector dlogn_bet = dlogn + bet.derive_cov( mp.flat_met_spher() ) ;
608 
609  source = qpig * a_car *( ener_euler + s_euler ) + ak_car
610  - dlogn(1)*dlogn_bet(1) - dlogn(2)*dlogn_bet(2) - dlogn(3)*dlogn_bet(3) ;
611  }
612  else {
613  source = qpig * nbar ;
614  }
615 
616  source.std_spectral_base() ;
617 
618  // Multiplication by -r^2 P_2(cos(theta))
619  // [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
620  // ------------------------------------------------------------------
621 
622  // Multiplication by r^2 :
623  // ----------------------
624  source.mult_r() ;
625  source.mult_r() ;
626  if (source.check_dzpuis(2)) {
627  source.inc_dzpuis(2) ;
628  }
629 
630  // Muliplication by cos^2(theta) :
631  // -----------------------------
632  Scalar temp = source ;
633 
634  // What follows is valid only for a mapping of class Map_radial :
635  assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ;
636 
637  if (temp.get_etat() == ETATQCQ) {
638  Valeur& vtemp = temp.set_spectral_va() ;
639  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
640  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
641  }
642 
643  // Muliplication by -P_2(cos(theta)) :
644  // ----------------------------------
645  source = 0.5 * source - 1.5 * temp ;
646 
647  // First part
648  // ------------
649 
650  p_mom_quad = new double( -source.integrale() / qpig ) ;
651 
652  // Correction terms from (3.37) of Friedman & Stergioulas' book
653  if (relativistic) {
654 
655  Scalar correct(mp) ;
656  correct = a_car * bbb * nn * press ;
657  correct.mult_rsint() ;
658  correct.std_spectral_base() ;
659 
660  double Bo = -32. * correct.integrale() / qpig ;
661 
662  *p_mom_quad -= 4./3. * (1./4. + Bo/(mass_g() * mass_g()))
663  * pow(mass_g(), 3) * ggrav * ggrav ;
664  }
665 
666  }
667 
668  return *p_mom_quad ;
669 
670 }
671 
672 
673  //----------------------------//
674  // Surface gravity //
675  //----------------------------//
676 
677 const Tbl& Star_rot::surf_grav() const {
678 
679  using namespace Unites ;
680  if (p_surf_grav == 0x0) { // a new computation is required
681 
682  const Mg3d* mg = mp.get_mg() ;
683  assert(mg->get_type_t() == SYM) ;
684  int nt = mg->get_nt(0) ; // number of theta points
685 
686  Scalar loggam = log(gam_euler) ;
687  loggam.std_spectral_base() ;
688 
690  - loggam.derive_cov( mp.flat_met_spher() );
691 
692  Scalar g_s = sqrt(contract( a, 0, a.up_down(gamma), 0 )) ;
693  g_s.std_spectral_base() ;
694 
695  p_surf_grav = new Tbl(nt) ;
697 
698  for(int j = 0; j < nt; j++) {
699  int l_b = l_surf()(0, j) ;
700  double xi_b = xi_surf()(0, j) ;
701  p_surf_grav->set(j) = g_s.get_spectral_va().val_point_jk(l_b, xi_b, j, 0) ;
702  }
703  }
704 
705  return *p_surf_grav ;
706 
707 }
708 
709  //----------------------------//
710  // Circular radius (\theta) //
711  //----------------------------//
712 
714 
715  using namespace Unites ;
716  if (p_radius_Morsink == 0x0) { // a new computation is required
717 
718  const Mg3d* mg = mp.get_mg() ;
719  assert(mg->get_type_t() == SYM) ;
720  int nt = mg->get_nt(0) ; // number of theta points
721 
722  p_radius_Morsink = new Tbl(nt) ;
724  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
725 
726  for(int j = 0; j < nt; j++) {
727  int l_b = l_surf()(0, j) ;
728  double xi_b = xi_surf()(0, j) ;
730  = bbb.get_spectral_va().val_point_jk(l_b, xi_b, j, 0)
731  * mp_rad->val_r_jk(l_b, xi_b, j, 0) ;
732  }
733  }
734 
735  return *p_radius_Morsink ;
736 
737 }
738 
739 }
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
double * p_mass_b
Baryon mass.
Definition: star.h:268
const Valeur & dsdt() const
Returns of *this.
Definition: valeur_dsdt.C:115
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
Scalar dzeta
Metric potential .
Definition: star_rot.h:149
virtual double mean_radius() const
Mean star radius from the area .
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:119
void mult_r()
Multiplication by r everywhere; dzpuis is not changed.
Map & mp
Mapping associated with the star.
Definition: star.h:180
static double lambda_grv2(const Scalar &sou_m, const Scalar &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
virtual double r_circ_merid() const
Crcumferential meridional radius.
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
double * p_r_circ
Circumferential radius (equator)
Definition: star_rot.h:251
Metric gamma
3-metric
Definition: star.h:235
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
Scalar bbb
Metric factor B.
Definition: star_rot.h:122
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
void annule_hard()
Sets the Valeur to zero in a hard way.
Definition: valeur.C:726
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
virtual double z_pole() const
Redshift factor at North pole.
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
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.
double integrale() const
Computes the integral over all space of *this .
Definition: scalar_integ.C:67
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition: grilles.h:504
Basic integer array class.
Definition: itbl.h:122
const Valeur & sx() const
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) ...
Definition: valeur_sx.C:113
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot.h:116
Scalar nbar
Baryon density in the fluid frame.
Definition: star.h:192
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:566
Scalar ent
Log-enthalpy.
Definition: star.h:190
virtual const Tbl & surf_grav() const
Surface gravity (table along the theta direction)
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition: star.h:201
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
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition: valeur.C:827
virtual double area() const
Integrated surface area in .
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...
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
const Tbl & xi_surf() const
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate ...
Definition: star_global.C:92
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:649
Tbl * p_xi_surf
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the su...
Definition: star.h:266
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
double * p_tsw
Ratio T/W.
Definition: star_rot.h:248
Scalar nphi
Metric coefficient .
Definition: star_rot.h:128
double * p_z_eqf
Forward redshift factor at equator.
Definition: star_rot.h:255
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: star.h:204
Scalar ener
Total energy density in the fluid frame.
Definition: star.h:193
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 void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Scalar press
Fluid pressure.
Definition: star.h:194
virtual double mass_b() const
Baryon mass.
const Valeur & ssint() const
Returns of *this.
Definition: valeur_ssint.C:115
const Mg3d * get_angu() const
Returns the pointer on the associated angular grid.
Definition: mg3d.C:609
virtual const Tbl & radius_Morsink() const
Physical radius following Morsink et al.
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:569
virtual double mom_quad() const
Quadrupole moment.
virtual double angu_mom() const
Angular momentum.
double * p_z_pole
Redshift factor at North pole.
Definition: star_rot.h:257
Base class for pure radial mappings.
Definition: map.h:1606
const Valeur & mult_st() const
Returns applied to *this.
virtual double val_r_jk(int l, double xi, int j, int k) const =0
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
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
double val_point_jk(int l, double x, int j, int k) const
Computes the value of the field represented by *this at an arbitrary point in , but collocation point...
Definition: valeur.C:903
Coord xsr
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain...
Definition: map.h:1619
void mult_rsint()
Multiplication by everywhere; dzpuis is not changed.
virtual double r_circ() const
Circumferential equatorial radius.
double * p_angu_mom
Angular momentum.
Definition: star_rot.h:247
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Scalar logn
Logarithm of the lapse N .
Definition: star.h:222
Scalar b_car
Square of the metric factor B.
Definition: star_rot.h:125
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition: scalar.C:896
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition: cmp.C:647
double ray_pole() const
Coordinate radius at [r_unit].
Definition: star_global.C:281
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
double * p_aplat
Flatening r_pole/r_eq.
Definition: star_rot.h:253
Multi-domain grid.
Definition: grilles.h:276
virtual double z_eqb() const
Backward redshift factor at equator.
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
const Scalar & srdsdt() const
Returns of *this .
Definition: scalar_deriv.C:144
Scalar nn
Lapse function N .
Definition: star.h:225
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:136
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:112
const Valeur & mult_ct() const
Returns applied to *this.
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:476
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:616
bool check_dzpuis(int dzi) const
Returns false if the last domain is compactified and *this is not zero in this domain and dzpuis is n...
Definition: scalar.C:879
Scalar ak_car
Scalar .
Definition: star_rot.h:201
virtual double aplat() const
Flatening r_pole/r_eq.
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:389
Tbl * p_surf_grav
Surface gravity (along the theta direction)
Definition: star_rot.h:267
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition: star.h:198
Itbl * p_l_surf
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surfa...
Definition: star.h:260
double * p_mass_g
Gravitational mass.
Definition: star.h:269
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:327
virtual double z_eqf() const
Forward redshift factor at equator.
virtual double tsw() const
Ratio T/W.
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:613
double * p_z_eqb
Backward redshift factor at equator.
Definition: star_rot.h:256