LORENE
et_rot_global.C
1 /*
2  * Methods for computing global quantities within the class Etoile_rot
3  *
4  * (see file etoile.h for documentation)
5  */
6 
7 /*
8  * Copyright (c) 2000-2001 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: et_rot_global.C,v 1.13 2026/02/26 09:53:04 j_novak Exp $
33  * $Log: et_rot_global.C,v $
34  * Revision 1.13 2026/02/26 09:53:04 j_novak
35  * Added method to compute and store circumferential radius as function of latitude (radius_Morsink).
36  *
37  * Revision 1.12 2023/07/04 08:59:53 j_novak
38  * Added method "r_circ_merid()" to compute the circumferential meridional radius.
39  *
40  * Revision 1.11 2017/10/06 12:36:34 a_sourie
41  * Cleaning of tabulated 2-fluid EoS class + superfluid rotating star model.
42  *
43  * Revision 1.10 2016/12/05 16:17:54 j_novak
44  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
45  *
46  * Revision 1.9 2015/06/12 12:38:25 j_novak
47  * Implementation of the corrected formula for the quadrupole momentum.
48  *
49  * Revision 1.8 2015/06/10 14:37:44 a_sourie
50  * Corrected the formula for the quadrupole.
51  *
52  * Revision 1.7 2014/10/13 08:52:57 j_novak
53  * Lorene classes and functions now belong to the namespace Lorene.
54  *
55  * Revision 1.6 2014/10/06 15:13:09 j_novak
56  * Modified #include directives to use c++ syntax.
57  *
58  * Revision 1.5 2012/08/12 17:48:35 p_cerda
59  * Magnetstar: New classes for magnetstar. Allowing for non-equatorial symmetry in Etoile et al. Adding B_phi in Et_rot_mag.
60  *
61  * Revision 1.4 2004/03/25 10:29:06 j_novak
62  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
63  *
64  * Revision 1.3 2003/11/03 16:47:13 e_gourgoulhon
65  * Corrected error in grv2() in the Newtonian case.
66  *
67  * Revision 1.2 2002/04/05 09:09:37 j_novak
68  * The inversion of the EOS for 2-fluids polytrope has been modified.
69  * Some errors in the determination of the surface were corrected.
70  *
71  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
72  * LORENE
73  *
74  * Revision 1.5 2000/11/19 18:52:09 eric
75  * grv2() operationnelle.
76  *
77  * Revision 1.4 2000/10/12 15:34:55 eric
78  * Calcul de la masse grav, de GRV3 et du moment quadrupolaire.
79  *
80  * Revision 1.3 2000/08/31 11:25:58 eric
81  * *** empty log message ***
82  *
83  * Revision 1.2 2000/08/25 12:28:16 eric
84  * *** empty log message ***
85  *
86  * Revision 1.1 2000/07/20 15:32:56 eric
87  * Initial revision
88  *
89  *
90  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_global.C,v 1.13 2026/02/26 09:53:04 j_novak Exp $
91  *
92  */
93 
94 // Headers C
95 #include <cstdlib>
96 #include <cmath>
97 
98 // Headers Lorene
99 #include "etoile.h"
100 #include "unites.h"
101 
102  //--------------------------//
103  // Stellar surface //
104  //--------------------------//
105 
106 namespace Lorene {
107 const Itbl& Etoile_rot::l_surf() const {
108 
109  if (p_l_surf == 0x0) { // a new computation is required
110 
111  assert(p_xi_surf == 0x0) ; // consistency check
112 
113  int np = mp.get_mg()->get_np(0) ;
114  int nt = mp.get_mg()->get_nt(0) ;
115 
116  p_l_surf = new Itbl(np, nt) ;
117  p_xi_surf = new Tbl(np, nt) ;
118 
119  double ent0 = 0 ; // definition of the surface
120  double precis = 1.e-15 ;
121  int nitermax = 100 ;
122  int niter ;
123 
124  (ent().va).equipot(ent0, nzet, precis, nitermax, niter, *p_l_surf,
125  *p_xi_surf) ;
126 
127  }
128 
129  return *p_l_surf ;
130 
131 }
132 
133  //--------------------------//
134  // Baryon mass //
135  //--------------------------//
136 
137 double Etoile_rot::mass_b() const {
138 
139  if (p_mass_b == 0x0) { // a new computation is required
140 
141  if (relativistic) {
142 
143  Cmp dens = a_car() * bbb() * gam_euler() * nbar() ;
144 
145  dens.std_base_scal() ;
146 
147  p_mass_b = new double( dens.integrale() ) ;
148 
149 
150  }
151  else{ // Newtonian case
152  assert(nbar.get_etat() == ETATQCQ) ;
153 
154  p_mass_b = new double( nbar().integrale() ) ;
155 
156  }
157 
158  }
159 
160  return *p_mass_b ;
161 
162 }
163 
164 
165  //----------------------------//
166  // Gravitational mass //
167  //----------------------------//
168 
169 double Etoile_rot::mass_g() const {
170 
171  if (p_mass_g == 0x0) { // a new computation is required
172 
173  if (relativistic) {
174 
175  Tenseur source = nnn * (ener_euler + s_euler)
176  + 2 * bbb * (ener_euler + press)
177  * tnphi * uuu ;
178  source = a_car * bbb * source ;
179 
180  source.set_std_base() ;
181 
182  p_mass_g = new double( source().integrale() ) ;
183 
184 
185  }
186  else{ // Newtonian case
187  p_mass_g = new double( mass_b() ) ; // in the Newtonian case
188  // M_g = M_b
189  }
190  }
191 
192  return *p_mass_g ;
193 
194 }
195 
196  //----------------------------//
197  // Angular momentum //
198  //----------------------------//
199 
200 double Etoile_rot::angu_mom() const {
201 
202  if (p_angu_mom == 0x0) { // a new computation is required
203 
204  Cmp dens = uuu() ;
205 
206  dens.mult_r() ; // Multiplication by
207  dens.va = (dens.va).mult_st() ; // r sin(theta)
208 
209  if (relativistic) {
210  dens = a_car() * b_car() * (ener_euler() + press())
211  * dens ;
212  }
213  else { // Newtonian case
214  dens = nbar() * dens ;
215  }
216 
217  dens.std_base_scal() ;
218 
219  p_angu_mom = new double( dens.integrale() ) ;
220 
221  }
222 
223  return *p_angu_mom ;
224 
225 }
226 
227 
228  //----------------------------//
229  // T/W //
230  //----------------------------//
231 
232 double Etoile_rot::tsw() const {
233 
234  if (p_tsw == 0x0) { // a new computation is required
235 
236  double tcin = 0.5 * omega * angu_mom() ;
237 
238  if (relativistic) {
239 
240  Cmp dens = a_car() * bbb() * gam_euler() * ener() ;
241  dens.std_base_scal() ;
242  double mass_p = dens.integrale() ;
243 
244  p_tsw = new double( tcin / ( mass_p + tcin - mass_g() ) ) ;
245 
246  }
247  else { // Newtonian case
248  Cmp dens = 0.5 * nbar() * logn() ;
249  dens.std_base_scal() ;
250  double wgrav = dens.integrale() ;
251  p_tsw = new double( tcin / fabs(wgrav) ) ;
252  }
253 
254 
255  }
256 
257  return *p_tsw ;
258 
259 }
260 
261 
262  //----------------------------//
263  // GRV2 //
264  //----------------------------//
265 
266 double Etoile_rot::grv2() const {
267 
268  using namespace Unites ;
269  if (p_grv2 == 0x0) { // a new computation is required
270 
271  Tenseur sou_m(mp) ;
272  if (relativistic) {
273  sou_m = 2 * qpig * a_car * (press + (ener_euler+press)
274  * uuu*uuu ) ;
275  }
276  else {
277  sou_m = 2 * qpig * (press + nbar * uuu*uuu ) ;
278  }
279 
280  Tenseur sou_q = 1.5 * ak_car
282  logn.gradient_spher() ) ;
283 
284  p_grv2 = new double( double(1) - lambda_grv2(sou_m(), sou_q()) ) ;
285 
286  }
287 
288  return *p_grv2 ;
289 
290 }
291 
292 
293  //----------------------------//
294  // GRV3 //
295  //----------------------------//
296 
297 double Etoile_rot::grv3(ostream* ost) const {
298 
299  using namespace Unites ;
300 
301  if (p_grv3 == 0x0) { // a new computation is required
302 
303  Tenseur source(mp) ;
304 
305  // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
306  // ------------------ Class. Quantum Grav. 11, 443 (1994)]
307 
308  if (relativistic) {
309  Tenseur alpha = dzeta - logn ;
310  Tenseur beta = log( bbb ) ;
311  beta.set_std_base() ;
312 
313  source = 0.75 * ak_car
315  logn.gradient_spher() )
316  + 0.5 * flat_scalar_prod(alpha.gradient_spher(),
317  beta.gradient_spher() ) ;
318 
319  Cmp aa = alpha() - 0.5 * beta() ;
320  Cmp daadt = aa.srdsdt() ; // 1/r d/dth
321 
322  // What follows is valid only for a mapping of class Map_radial :
323  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
324  if (mpr == 0x0) {
325  cout << "Etoile_rot::grv3: the mapping does not belong"
326  << " to the class Map_radial !" << endl ;
327  abort() ;
328  }
329 
330  // Computation of 1/tan(theta) * 1/r daa/dtheta
331  if (daadt.get_etat() == ETATQCQ) {
332  Valeur& vdaadt = daadt.va ;
333  vdaadt = vdaadt.ssint() ; // division by sin(theta)
334  vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
335  }
336 
337  Cmp temp = aa.dsdr() + daadt ;
338  temp = ( bbb() - a_car()/bbb() ) * temp ;
339  temp.std_base_scal() ;
340 
341  // Division by r
342  Valeur& vtemp = temp.va ;
343  vtemp = vtemp.sx() ; // division by xi in the nucleus
344  // Id in the shells
345  // division by xi-1 in the ZEC
346  vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
347  // by 1/r in the shells
348  // by r(xi-1) in the ZEC
349 
350  // In the ZEC, a multiplication by r has been performed instead
351  // of the division:
352  temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;
353 
354  source = bbb() * source() + 0.5 * temp ;
355 
356  }
357  else{
358  source = - 0.5 * flat_scalar_prod(logn.gradient_spher(),
359  logn.gradient_spher() ) ;
360  }
361 
362  source.set_std_base() ;
363 
364  double int_grav = source().integrale() ;
365 
366  // Matter term
367  // -----------
368 
369  if (relativistic) {
370  source = qpig * a_car * bbb * s_euler ;
371  }
372  else{
373  source = qpig * ( 3 * press + nbar * uuu * uuu ) ;
374  }
375 
376  source.set_std_base() ;
377 
378  double int_mat = source().integrale() ;
379 
380  // Virial error
381  // ------------
382  if (ost != 0x0) {
383  *ost << "Etoile_rot::grv3 : gravitational term : " << int_grav
384  << endl ;
385  *ost << "Etoile_rot::grv3 : matter term : " << int_mat
386  << endl ;
387  }
388 
389  p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;
390 
391  }
392 
393  return *p_grv3 ;
394 
395 }
396 
397 
398  //----------------------------//
399  // R_circ //
400  //----------------------------//
401 
402 double Etoile_rot::r_circ() const {
403 
404  if (p_r_circ == 0x0) { // a new computation is required
405 
406  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
407  const Mg3d* mg = mp.get_mg() ;
408 
409 
410  int l_b = nzet - 1 ;
411  int i_b = mg->get_nr(l_b) - 1 ;
412  int j_b;
413  if (mg->get_type_t() == SYM) {
414  j_b = mg->get_nt(l_b) - 1 ;
415  }else{
416  j_b = (mg->get_nt(l_b) - 1)/2 ;
417  }
418  int k_b = 0 ;
419 
420  p_r_circ = new double( bbb()(l_b, k_b, j_b, i_b) * ray_eq() ) ;
421 
422  }
423 
424  return *p_r_circ ;
425 
426 }
427 
428 double Etoile_rot::r_circ_merid() const {
429 
430  if (p_r_circ_merid == 0x0) { // a new computation is required
431  const Mg3d& mg = *(mp.get_mg()) ;
432  int np = mg.get_np(0) ;
433  int nt = mg.get_nt(0) ;
434  assert(np == 1) ; //Only valid for axisymmetric configurations
435 
436  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
437  assert(mp_rad != 0x0) ;
438 
439  Valeur va_r(mg.get_angu()) ;
440  va_r.annule_hard() ;
441  Itbl lsurf = l_surf() ;
442  Tbl xisurf = xi_surf() ;
443  for (int k=0; k<np; k++) {
444  for (int j=0; j<nt; j++) {
445  int l_star = lsurf(k, j) ;
446  double xi_star = xisurf(k, j) ;
447  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star, xi_star, j, k) ;
448  }
449  }
450  va_r.std_base_scal() ;
451 
452  Valeur integ(mg.get_angu()) ;
453  Valeur dr = va_r.dsdt() ;
454  integ = sqrt(va_r*va_r + dr*dr) ;
455  Cmp aaaa = get_a_car()() ;
456  Valeur a2 = aaaa.va ; a2.std_base_scal() ;
457  int k = 0 ;
458  for (int j=0; j<nt; j++) {
459  integ.set(0, k, j, 0)
460  *= sqrt(a2.val_point_jk(lsurf(k, j), xisurf(k, j), j, k)) ;
461  }
462  integ.std_base_scal() ;
463  integ.coef() ;
464  double surftot = (*integ.c_cf)(0, k, 0, 0) ; // Only j=0 has a nonvanishing
465  // integral
466  p_r_circ_merid = new double(surftot ) ;
467 
468  }
469 
470  return *p_r_circ_merid ;
471 
472 }
473 
474 
475  //----------------------------//
476  // Surface area //
477  //----------------------------//
478 
479  double Etoile_rot::area() const {
480 
481  if (p_area == 0x0) { // a new computation is required
482  const Mg3d& mg = *(mp.get_mg()) ;
483  int np = mg.get_np(0) ;
484  int nt = mg.get_nt(0) ;
485  assert(np == 1) ; //Only valid for axisymmetric configurations
486 
487  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
488  assert(mp_rad != 0x0) ;
489 
490  Valeur va_r(mg.get_angu()) ;
491  va_r.annule_hard() ;
492  Itbl lsurf = l_surf() ;
493  Tbl xisurf = xi_surf() ;
494 
495  for (int k=0; k<np; k++) {
496  for (int j=0; j<nt; j++) {
497  int l_star = lsurf(k, j) ;
498  double xi_star = xisurf(k, j) ;
499 
500  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star, xi_star, j, k) ;
501  }
502  }
503  va_r.std_base_scal() ;
504 
505  Valeur integ(mg.get_angu()) ;
506  Valeur dr = va_r.dsdt() ;
507  integ = sqrt(va_r*va_r + dr*dr) ;
508  Cmp aaaa = get_a_car()() ;
509  Valeur a2 = aaaa.va ; a2.std_base_scal() ;
510  Cmp bbbb = get_bbb()() ;
511  Valeur b = bbbb.va ; b.std_base_scal() ;
512  for (int k=0; k<np; k++) {
513  for (int j=0; j<nt; j++) {
514  integ.set(0, k, j, 0) *= sqrt(a2.val_point_jk(lsurf(k, j), xisurf(k, j), j, k))
515  * b.val_point_jk(lsurf(k, j), xisurf(k, j), j, k) * va_r(0, k, j, 0) ;
516  }
517  }
518  integ.std_base_scal() ;
519  Valeur integ2 = integ.mult_st() ;
520  double surftot = 0. ;
521  for (int j=0; j<nt; j++) {
522  surftot += (*integ2.c_cf)(0, 0, j, 0) / double(2*j+1) ;
523  }
524 
525  p_area = new double( 4*M_PI*surftot) ;
526 
527  }
528 
529  return *p_area ;
530 
531 }
532 
533  double Etoile_rot::mean_radius() const {
534 
535  return sqrt(area()/(4*M_PI)) ;
536 
537  }
538 
539 
540 
541 
542 
543 
544 
545 
546 
547 
548 
549 
550 
551 
552 
553  //----------------------------//
554  // Flattening //
555  //----------------------------//
556 
557 double Etoile_rot::aplat() const {
558 
559  if (p_aplat == 0x0) { // a new computation is required
560 
561  p_aplat = new double( ray_pole() / ray_eq() ) ;
562 
563  }
564 
565  return *p_aplat ;
566 
567 }
568 
569 
570  //----------------------------//
571  // Z_eq_f //
572  //----------------------------//
573 
574 double Etoile_rot::z_eqf() const {
575 
576  if (p_z_eqf == 0x0) { // a new computation is required
577 
578  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
579  const Mg3d* mg = mp.get_mg() ;
580  int l_b = nzet - 1 ;
581  int i_b = mg->get_nr(l_b) - 1 ;
582  int j_b;
583  if (mg->get_type_t() == SYM) {
584  j_b = mg->get_nt(l_b) - 1 ;
585  }else{
586  j_b = (mg->get_nt(l_b) - 1)/2 ;
587  }
588  int k_b = 0 ;
589 
590  double u_eq = uuu()(l_b, k_b, j_b, i_b) ;
591  double n_eq = nnn()(l_b, k_b, j_b, i_b) ;
592  double nphi_eq = nphi()(l_b, k_b, j_b, i_b) ;
593 
594  p_z_eqf = new double( sqrt((1.-u_eq)/(1.+u_eq))
595  / (n_eq + nphi_eq * r_circ() )
596  - 1. ) ;
597  }
598 
599  return *p_z_eqf ;
600 
601 }
602  //----------------------------//
603  // Z_eq_b //
604  //----------------------------//
605 
606 double Etoile_rot::z_eqb() const {
607 
608  if (p_z_eqb == 0x0) { // a new computation is required
609 
610  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
611  const Mg3d* mg = mp.get_mg() ;
612  int l_b = nzet - 1 ;
613  int i_b = mg->get_nr(l_b) - 1 ;
614  int j_b;
615  if (mg->get_type_t() == SYM) {
616  j_b = mg->get_nt(l_b) - 1 ;
617  }else{
618  j_b = (mg->get_nt(l_b) - 1) / 2 ;
619  }
620  int k_b = 0 ;
621 
622  double u_eq = uuu()(l_b, k_b, j_b, i_b) ;
623  double n_eq = nnn()(l_b, k_b, j_b, i_b) ;
624  double nphi_eq = nphi()(l_b, k_b, j_b, i_b) ;
625 
626  p_z_eqb = new double( sqrt((1.+u_eq)/(1.-u_eq))
627  / (n_eq - nphi_eq * r_circ() )
628  - 1. ) ;
629 
630  }
631 
632  return *p_z_eqb ;
633 
634 }
635 
636 
637  //----------------------------//
638  // Z_pole //
639  //----------------------------//
640 
641 double Etoile_rot::z_pole() const {
642 
643  if (p_z_pole == 0x0) { // a new computation is required
644 
645  double n_pole = nnn().val_point(ray_pole(), 0., 0.) ;
646 
647  p_z_pole = new double( 1. / n_pole - 1. ) ;
648 
649  }
650 
651  return *p_z_pole ;
652 
653 }
654 
655 
656  //----------------------------//
657  // Quadrupole moment //
658  //----------------------------//
659 
660 
661 double Etoile_rot::mom_quad() const {
662 
663  using namespace Unites ;
664 
665  if (p_mom_quad == 0x0) { // a new computation is required
666 
667  p_mom_quad = new double( mom_quad_old() ) ;
668  if (relativistic) {
669  double b = mom_quad_Bo() / ( mass_g() * mass_g() ) ;
670  *p_mom_quad -= 4./3. * ( 1./4. + b ) * pow(mass_g(), 3) * ggrav * ggrav ;
671  }
672  }
673 
674  return *p_mom_quad ;
675 
676 }
677 
678 
679 double Etoile_rot::mom_quad_Bo() const {
680 
681  using namespace Unites ;
682 
683  if (p_mom_quad_Bo == 0x0) { // a new computation is required
684 
685  Cmp dens(mp) ;
686 
687  dens = press() ;
688  dens = a_car() * bbb() * nnn() * dens ;
689  dens.mult_rsint() ;
690  dens.std_base_scal() ;
691 
692  p_mom_quad_Bo = new double( - 32. * dens.integrale() / qpig ) ;
693 
694  }
695 
696  return *p_mom_quad_Bo ;
697 
698 }
699 
700 
701 
702 double Etoile_rot::mom_quad_old() const {
703 
704  using namespace Unites ;
705 
706  if (p_mom_quad_old == 0x0) { // a new computation is required
707 
708  // Source for of the Poisson equation for nu
709  // -----------------------------------------
710 
711  Tenseur source(mp) ;
712 
713  if (relativistic) {
714  Tenseur beta = log(bbb) ;
715  beta.set_std_base() ;
716  source = qpig * a_car *( ener_euler + s_euler )
718  logn.gradient_spher() + beta.gradient_spher()) ;
719  }
720  else {
721  source = qpig * nbar ;
722  }
723  source.set_std_base() ;
724 
725  // Multiplication by -r^2 P_2(cos(theta))
726  // [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
727  // ------------------------------------------------------------------
728 
729  // Multiplication by r^2 :
730  // ----------------------
731  Cmp& csource = source.set() ;
732  csource.mult_r() ;
733  csource.mult_r() ;
734  if (csource.check_dzpuis(2)) {
735  csource.inc2_dzpuis() ;
736  }
737 
738  // Muliplication by cos^2(theta) :
739  // -----------------------------
740  Cmp temp = csource ;
741 
742  // What follows is valid only for a mapping of class Map_radial :
743  assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ;
744 
745  if (temp.get_etat() == ETATQCQ) {
746  Valeur& vtemp = temp.va ;
747  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
748  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
749  }
750 
751  // Muliplication by -P_2(cos(theta)) :
752  // ----------------------------------
753  source = 0.5 * source() - 1.5 * temp ;
754 
755  // Final result
756  // ------------
757 
758  p_mom_quad_old = new double(- source().integrale() / qpig ) ;
759 
760  }
761 
762  return *p_mom_quad_old ;
763 
764 }
765 
766  //----------------------------//
767  // Circular radius (\theta) //
768  //----------------------------//
769 
771 
772  using namespace Unites ;
773  if (p_radius_Morsink == 0x0) { // a new computation is required
774 
775  const Mg3d* mg = mp.get_mg() ;
776  assert(mg->get_type_t() == SYM) ;
777  int nt = mg->get_nt(0) ; // number of theta points
778 
779  p_radius_Morsink = new Tbl(nt) ;
781  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
782 
783  for(int j = 0; j < nt; j++) {
784  int l_b = l_surf()(0, j) ;
785  double xi_b = xi_surf()(0, j) ;
787  = bbb().va.val_point_jk(l_b, xi_b, j, 0)
788  * mp_rad->val_r_jk(l_b, xi_b, j, 0) ;
789  }
790  }
791 
792  return *p_radius_Morsink ;
793 
794 }
795 
796 
797 
798 
799 }
double * p_z_eqb
Backward redshift factor at equator.
Definition: etoile.h:1646
virtual double z_eqf() const
Forward redshift factor at equator.
const Cmp & dsdr() const
Returns of *this .
Definition: cmp_deriv.C:87
double * p_z_pole
Redshift factor at North pole.
Definition: etoile.h:1647
virtual double mom_quad_Bo() const
Part of the quadrupole moment.
virtual double r_circ() const
Circumferential equatorial radius.
const Tbl & xi_surf() const
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate...
double * p_mom_quad_old
Part of the quadrupole moment.
Definition: etoile.h:1649
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
Tbl * p_radius_Morsink
Physical radius close to the definition of Morsink et al. (2007)
Definition: etoile.h:1659
double * p_r_circ
Circumferential equatorial radius.
Definition: etoile.h:1641
double * p_mom_quad_Bo
Part of the quadrupole moment.
Definition: etoile.h:1650
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:449
Tenseur tnphi
Component of the shift vector.
Definition: etoile.h:1521
const Tenseur & get_a_car() const
Returns the total conformal factor .
Definition: etoile.h:739
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition: tenseur.C:1564
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
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
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:791
int get_etat() const
Returns the logical state.
Definition: cmp.h:902
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tenseur nnn
Total lapse function.
Definition: etoile.h:515
double * p_z_eqf
Forward redshift factor at equator.
Definition: etoile.h:1645
Tenseur nphi
Metric coefficient .
Definition: etoile.h:1516
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:474
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
double * p_mom_quad
Quadrupole moment.
Definition: etoile.h:1648
const Cmp & srdsdt() const
Returns of *this .
Definition: cmp_deriv.C:108
Tenseur b_car
Square of the metric factor B.
Definition: etoile.h:1513
virtual double mom_quad_old() const
Part of the quadrupole moment.
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
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
Tenseur press
Fluid pressure.
Definition: etoile.h:467
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...
virtual double mean_radius() const
Mean radius.
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition: valeur.C:827
double * p_aplat
Flatening r_pole/r_eq.
Definition: etoile.h:1644
virtual double r_circ_merid() const
Circumferential meridional radius.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
void mult_r()
Multiplication by r everywhere.
Definition: cmp_r_manip.C:94
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
double * p_r_circ_merid
Circumferential meridional radius.
Definition: etoile.h:1642
double * p_grv3
Error on the virial identity GRV3.
Definition: etoile.h:1640
virtual double angu_mom() const
Angular momentum.
const Valeur & ssint() const
Returns of *this.
Definition: valeur_ssint.C:115
virtual double z_pole() const
Redshift factor at North pole.
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:465
const Mg3d * get_angu() const
Returns the pointer on the associated angular grid.
Definition: mg3d.C:609
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:477
double * p_tsw
Ratio T/W.
Definition: etoile.h:1638
virtual double grv2() const
Error on the virial identity GRV2.
virtual double mass_b() const
Baryon mass.
double * p_area
Surface area.
Definition: etoile.h:1643
Base class for pure radial mappings.
Definition: map.h:1611
void mult_rsint()
Multiplication by .
Definition: cmp_r_manip.C:119
Map & mp
Mapping associated with the star.
Definition: etoile.h:435
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...
double integrale() const
Computes the integral over all space of *this .
Definition: cmp_integ.C:58
Tbl * p_xi_surf
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the su...
Definition: etoile.h:551
virtual double z_eqb() const
Backward redshift factor at equator.
Tenseur ak_car
Scalar .
Definition: etoile.h:1592
int get_etat() const
Returns the logical state.
Definition: tenseur.h:713
static double lambda_grv2(const Cmp &sou_m, const Cmp &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
virtual double tsw() const
Ratio T/W.
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:1624
Tenseur bbb
Metric factor B.
Definition: etoile.h:1510
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
void inc2_dzpuis()
Increases by 2 the value of dzpuis and changes accordingly the values of the Cmp in the external comp...
Definition: cmp_r_manip.C:195
Tenseur uuu
Norm of u_euler.
Definition: etoile.h:1524
double * p_mass_g
Gravitational mass.
Definition: etoile.h:554
double omega
Rotation angular velocity ([f_unit] )
Definition: etoile.h:1507
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition: cmp.C:647
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:438
virtual double mom_quad() const
Quadrupole moment.
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
Tenseur a_car
Total conformal factor .
Definition: etoile.h:521
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:443
Multi-domain grid.
Definition: grilles.h:276
double ray_pole() const
Coordinate radius at [r_unit].
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:466
int get_dzpuis() const
Returns dzpuis.
Definition: cmp.h:906
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:463
virtual double area() const
Surface area.
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: cmp.C:718
virtual const Tbl & radius_Morsink() const
Physical radius following Morsink et al.
const Valeur & mult_ct() const
Returns applied to *this.
Tenseur & logn
Metric potential = logn_auto.
Definition: etoile.h:1527
void set_dzpuis(int)
Set a value to dzpuis.
Definition: cmp.C:657
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
virtual double mass_g() const
Gravitational mass.
double * p_grv2
Error on the virial identity GRV2.
Definition: etoile.h:1639
const Tenseur & get_bbb() const
Returns the metric factor B.
Definition: etoile.h:1722
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:471
Tenseur & dzeta
Metric potential = beta_auto.
Definition: etoile.h:1540
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:467
double * p_mass_b
Baryon mass.
Definition: etoile.h:553
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:307
double * p_angu_mom
Angular momentum.
Definition: etoile.h:1637
virtual double aplat() const
Flatening r_pole/r_eq.
Itbl * p_l_surf
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surfa...
Definition: etoile.h:545