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.11 2023/07/04 09:01:09 j_novak Exp $
33  * $Log: star_rot_global.C,v $
34  * Revision 1.11 2023/07/04 09:01:09 j_novak
35  * Changed the name r_pol_circ -> r_circ_merid to explicitly show the radius is a circumferential one, along a meridian.
36  *
37  * Revision 1.10 2023/07/03 14:32:31 j_novak
38  * Added method "r_pol_circ()", to compute polar circumferential radius, with a circumference along a meridian.
39  *
40  * Revision 1.9 2021/10/06 09:41:28 j_novak
41  * Corrected quadrupole computation from Friedman & Stergioulas book
42  *
43  * Revision 1.8 2017/04/13 09:59:34 m_bejger
44  * Bug fix (contraction using physical gamma metric in surf_grav - vector a)
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 2016/12/05 16:18:15 j_novak
50  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
51  *
52  * Revision 1.5 2015/05/19 09:30:56 j_novak
53  * New methods for computing the area of the star and its mean radius.
54  *
55  * Revision 1.4 2014/10/13 08:53:39 j_novak
56  * Lorene classes and functions now belong to the namespace Lorene.
57  *
58  * Revision 1.3 2014/10/06 15:13:17 j_novak
59  * Modified #include directives to use c++ syntax.
60  *
61  * Revision 1.2 2010/01/25 22:33:35 e_gourgoulhon
62  * Debugging...
63  *
64  * Revision 1.1 2010/01/25 18:15:52 e_gourgoulhon
65  * First version.
66  *
67  *
68  *
69  * $Header: /cvsroot/Lorene/C++/Source/Star/star_rot_global.C,v 1.11 2023/07/04 09:01:09 j_novak Exp $
70  *
71  */
72 
73 // Headers C
74 #include <cstdlib>
75 #include <cmath>
76 
77 // Headers Lorene
78 #include "star_rot.h"
79 #include "unites.h"
80 
81  //--------------------------//
82  // Stellar surface //
83  //--------------------------//
84 
85 namespace Lorene {
86  const Itbl& Star_rot::l_surf() const {
87 
88  if (p_l_surf == 0x0) { // a new computation is required
89 
90  assert(p_xi_surf == 0x0) ; // consistency check
91 
92  int np = mp.get_mg()->get_np(0) ;
93  int nt = mp.get_mg()->get_nt(0) ;
94 
95  p_l_surf = new Itbl(np, nt) ;
96  p_xi_surf = new Tbl(np, nt) ;
97 
98  double ent0 = 0 ; // definition of the surface
99  double precis = 1.e-15 ;
100  int nitermax = 100 ;
101  int niter ;
102 
103  (ent.get_spectral_va()).equipot(ent0, nzet, precis, nitermax, niter, *p_l_surf,
104  *p_xi_surf) ;
105 
106  }
107 
108  return *p_l_surf ;
109 
110  }
111 
112  //------------------------------//
113  // Baryon mass //
114  //------------------------------//
115 
116  double Star_rot::mass_b() const {
117 
118  if (p_mass_b == 0x0) { // a new computation is required
119 
120  if (relativistic) {
121 
122  Scalar dens = a_car * bbb * gam_euler * nbar ;
123 
124  p_mass_b = new double( dens.integrale() ) ;
125  }
126  else{ // Newtonian case
127  assert(nbar.get_etat() == ETATQCQ) ;
128 
129  p_mass_b = new double( nbar.integrale() ) ;
130 
131  }
132 
133  }
134  return *p_mass_b ;
135  }
136 
137 
138  //----------------------------//
139  // Gravitational mass //
140  //----------------------------//
141 
142  double Star_rot::mass_g() const {
143 
144  if (p_mass_g == 0x0) { // a new computation is required
145 
146  if (relativistic) {
147 
148  Scalar source = nn * (ener_euler + s_euler)
149  + 2 * bbb * (ener_euler + press)
150  * tnphi * uuu ;
151  source = a_car * bbb * source ;
152  source.std_spectral_base() ;
153 
154  p_mass_g = new double( source.integrale() ) ;
155  }
156  else{ // Newtonian case
157  p_mass_g = new double( mass_b() ) ; // in the Newtonian case
158  // M_g = M_b
159  }
160  }
161  return *p_mass_g ;
162 }
163 
164  //----------------------------//
165  // Angular momentum //
166  //----------------------------//
167 
168  double Star_rot::angu_mom() const {
169 
170  if (p_angu_mom == 0x0) { // a new computation is required
171 
172  Scalar dens = uuu ;
173 
174  dens.mult_r() ; // Multiplication by
175  dens.set_spectral_va() = (dens.get_spectral_va()).mult_st() ; // r sin(theta)
176 
177  if (relativistic) {
178  dens = a_car * b_car * (ener_euler + press) * dens ;
179  }
180  else { // Newtonian case
181  dens = nbar * dens ;
182  }
183 
184  p_angu_mom = new double( dens.integrale() ) ;
185  }
186  return *p_angu_mom ;
187  }
188 
189 
190  //----------------------------//
191  // T/W //
192  //----------------------------//
193 
194  double Star_rot::tsw() const {
195 
196  if (p_tsw == 0x0) { // a new computation is required
197 
198  double tcin = 0.5 * omega * angu_mom() ;
199 
200  if (relativistic) {
201 
202  Scalar dens = a_car * bbb * gam_euler * ener ;
203  double mass_p = dens.integrale() ;
204 
205  p_tsw = new double( tcin / ( mass_p + tcin - mass_g() ) ) ;
206 
207  }
208  else { // Newtonian case
209  Scalar dens = 0.5 * nbar * logn ;
210  double wgrav = dens.integrale() ;
211  p_tsw = new double( tcin / fabs(wgrav) ) ;
212  }
213  }
214  return *p_tsw ;
215 }
216 
217 
218  //----------------------------//
219  // GRV2 //
220  //----------------------------//
221 
222 double Star_rot::grv2() const {
223 
224  using namespace Unites ;
225  if (p_grv2 == 0x0) { // a new computation is required
226 
227  Scalar sou_m(mp) ;
228  if (relativistic) {
229  sou_m = 2 * qpig * a_car * (press + (ener_euler+press)
230  * uuu*uuu ) ;
231  }
232  else {
233  sou_m = 2 * qpig * (press + nbar * uuu*uuu ) ;
234  }
235 
236  Vector dlogn = logn.derive_cov( mp.flat_met_spher() ) ;
237  Scalar sou_q = 1.5 * ak_car
238  - dlogn(1)*dlogn(1) - dlogn(2)*dlogn(2) - dlogn(3)*dlogn(3) ;
239 
240  p_grv2 = new double( double(1) - lambda_grv2(sou_m, sou_q) ) ;
241 
242  }
243 
244  return *p_grv2 ;
245 
246 }
247 
248 
249  //----------------------------//
250  // GRV3 //
251  //----------------------------//
252 
253 double Star_rot::grv3(ostream* ost) const {
254 
255  using namespace Unites ;
256 
257  if (p_grv3 == 0x0) { // a new computation is required
258 
259  Scalar source(mp) ;
260 
261  // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
262  // ------------------ Class. Quantum Grav. 11, 443 (1994)]
263 
264  Vector dlogn = logn.derive_cov( mp.flat_met_spher() ) ;
265 
266  if (relativistic) {
267  Scalar alpha = dzeta - logn ;
268  Scalar bet = log( bbb ) ;
269  bet.std_spectral_base() ;
270 
271  Vector dalpha = alpha.derive_cov( mp.flat_met_spher() ) ;
272  Vector dbet = bet.derive_cov( mp.flat_met_spher() ) ;
273 
274  source = 0.75 * ak_car
275  - dlogn(1)*dlogn(1) - dlogn(2)*dlogn(2) - dlogn(3)*dlogn(3)
276  + 0.5 * ( dalpha(1)*dbet(1) + dalpha(2)*dbet(2) + dalpha(3)*dbet(3) ) ;
277 
278  Scalar aa = alpha - 0.5 * bet ;
279  Scalar daadt = aa.srdsdt() ; // 1/r d/dth
280 
281  // What follows is valid only for a mapping of class Map_radial :
282  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
283  if (mpr == 0x0) {
284  cout << "Star_rot::grv3: the mapping does not belong"
285  << " to the class Map_radial !" << endl ;
286  abort() ;
287  }
288 
289  // Computation of 1/tan(theta) * 1/r daa/dtheta
290  if (daadt.get_etat() == ETATQCQ) {
291  Valeur& vdaadt = daadt.set_spectral_va() ;
292  vdaadt = vdaadt.ssint() ; // division by sin(theta)
293  vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
294  }
295 
296  Scalar temp = aa.dsdr() + daadt ;
297  temp = ( bbb - a_car/bbb ) * temp ;
298  temp.std_spectral_base() ;
299 
300  // Division by r
301  Valeur& vtemp = temp.set_spectral_va() ;
302  vtemp = vtemp.sx() ; // division by xi in the nucleus
303  // Id in the shells
304  // division by xi-1 in the ZEC
305  vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
306  // by 1/r in the shells
307  // by r(xi-1) in the ZEC
308 
309  // In the ZEC, a multiplication by r has been performed instead
310  // of the division:
311  temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;
312 
313  source = bbb * source + 0.5 * temp ;
314 
315  }
316  else{
317  source = - 0.5 * ( dlogn(1)*dlogn(1) + dlogn(2)*dlogn(2) + dlogn(3)*dlogn(3) ) ;
318  }
319 
320  source.std_spectral_base() ;
321 
322  double int_grav = source.integrale() ;
323 
324  // Matter term
325  // -----------
326 
327  if (relativistic) {
328  source = qpig * a_car * bbb * s_euler ;
329  }
330  else{
331  source = qpig * ( 3 * press + nbar * uuu * uuu ) ;
332  }
333 
334  source.std_spectral_base() ;
335 
336  double int_mat = source.integrale() ;
337 
338  // Virial error
339  // ------------
340  if (ost != 0x0) {
341  *ost << "Star_rot::grv3 : gravitational term : " << int_grav
342  << endl ;
343  *ost << "Star_rot::grv3 : matter term : " << int_mat
344  << endl ;
345  }
346 
347  p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;
348 
349  }
350 
351  return *p_grv3 ;
352 
353 }
354 
355 
356  //----------------------------//
357  // R_circ //
358  //----------------------------//
359 
360 double Star_rot::r_circ() const {
361 
362  if (p_r_circ == 0x0) { // a new computation is required
363 
364  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
365  const Mg3d* mg = mp.get_mg() ;
366  assert(mg->get_type_t() == SYM) ;
367  int l_b = nzet - 1 ;
368  int i_b = mg->get_nr(l_b) - 1 ;
369  int j_b = mg->get_nt(l_b) - 1 ;
370  int k_b = 0 ;
371 
372  p_r_circ = new double( bbb.val_grid_point(l_b, k_b, j_b, i_b) * ray_eq() ) ;
373 
374  }
375 
376  return *p_r_circ ;
377 
378 }
379 
380 double Star_rot::r_circ_merid() const {
381 
382  if (p_r_circ_merid == 0x0) { // a new computation is required
383  const Mg3d& mg = *(mp.get_mg()) ;
384  int np = mg.get_np(0) ;
385  int nt = mg.get_nt(0) ;
386  assert(np == 1) ; //Only valid for axisymmetric configurations
387 
388  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
389  assert(mp_rad != 0x0) ;
390 
391  Valeur va_r(mg.get_angu()) ;
392  va_r.annule_hard() ;
393  Itbl lsurf = l_surf() ;
394  Tbl xisurf = xi_surf() ;
395  for (int k=0; k<np; k++) {
396  for (int j=0; j<nt; j++) {
397  int l_star = lsurf(k, j) ;
398  double xi_star = xisurf(k, j) ;
399  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star, xi_star, j, k) ;
400  }
401  }
402  va_r.std_base_scal() ;
403 
404  Valeur integ(mg.get_angu()) ;
405  Valeur dr = va_r.dsdt() ;
406  integ = sqrt(va_r*va_r + dr*dr) ;
408  int k = 0 ;
409  for (int j=0; j<nt; j++) {
410  integ.set(0, k, j, 0)
411  *= sqrt(a2.val_point_jk(lsurf(k, j), xisurf(k, j), j, k)) ;
412  }
413  integ.std_base_scal() ;
414  integ.coef() ;
415  double surftot = (*integ.c_cf)(0, k, 0, 0) ; // Only j=0 has a nonvanishing
416  // integral
417  p_r_circ_merid = new double(surftot ) ;
418 
419  }
420 
421  return *p_r_circ_merid ;
422 
423 }
424 
425  //----------------------------//
426  // Surface area //
427  //----------------------------//
428 
429  double Star_rot::area() const {
430 
431  if (p_area == 0x0) { // a new computation is required
432  const Mg3d& mg = *(mp.get_mg()) ;
433  int np = mg.get_np(0) ;
434  int nt = mg.get_nt(0) ;
435  assert(np == 1) ; //Only valid for axisymmetric configurations
436 
437  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
438  assert(mp_rad != 0x0) ;
439 
440  Valeur va_r(mg.get_angu()) ;
441  va_r.annule_hard() ;
442  Itbl lsurf = l_surf() ;
443  Tbl xisurf = xi_surf() ;
444  for (int k=0; k<np; k++) {
445  for (int j=0; j<nt; j++) {
446  int l_star = lsurf(k, j) ;
447  double xi_star = xisurf(k, j) ;
448  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star, xi_star, j, k) ;
449  }
450  }
451  va_r.std_base_scal() ;
452 
453  Valeur integ(mg.get_angu()) ;
454  Valeur dr = va_r.dsdt() ;
455  integ = sqrt(va_r*va_r + dr*dr) ;
458  for (int k=0; k<np; k++) {
459  for (int j=0; j<nt; j++) {
460  integ.set(0, k, j, 0) *= sqrt(a2.val_point_jk(lsurf(k, j), xisurf(k, j), j, k))
461  * b.val_point_jk(lsurf(k, j), xisurf(k, j), j, k) * va_r(0, k, j, 0) ;
462  }
463  }
464  integ.std_base_scal() ;
465  Valeur integ2 = integ.mult_st() ;
466  double surftot = 0. ;
467  for (int j=0; j<nt; j++) {
468  surftot += (*integ2.c_cf)(0, 0, j, 0) / double(2*j+1) ;
469  }
470 
471  p_area = new double( 4*M_PI*surftot) ;
472 
473  }
474 
475  return *p_area ;
476 
477 }
478 
479  double Star_rot::mean_radius() const {
480 
481  return sqrt(area()/(4*M_PI)) ;
482 
483  }
484 
485  //----------------------------//
486  // Flattening //
487  //----------------------------//
488 
489 double Star_rot::aplat() const {
490 
491  if (p_aplat == 0x0) { // a new computation is required
492 
493  p_aplat = new double( ray_pole() / ray_eq() ) ;
494 
495  }
496 
497  return *p_aplat ;
498 
499 }
500 
501 
502  //----------------------------//
503  // Z_eq_f //
504  //----------------------------//
505 
506 double Star_rot::z_eqf() const {
507 
508  if (p_z_eqf == 0x0) { // a new computation is required
509 
510  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
511  const Mg3d* mg = mp.get_mg() ;
512  assert(mg->get_type_t() == SYM) ;
513  int l_b = nzet - 1 ;
514  int i_b = mg->get_nr(l_b) - 1 ;
515  int j_b = mg->get_nt(l_b) - 1 ;
516  int k_b = 0 ;
517 
518  double u_eq = uuu.val_grid_point(l_b, k_b, j_b, i_b) ;
519  double n_eq = nn.val_grid_point(l_b, k_b, j_b, i_b) ;
520  double nphi_eq = nphi.val_grid_point(l_b, k_b, j_b, i_b) ;
521 
522  p_z_eqf = new double( sqrt((1.-u_eq)/(1.+u_eq))
523  / (n_eq + nphi_eq * r_circ() )
524  - 1. ) ;
525  }
526 
527  return *p_z_eqf ;
528 
529 }
530  //----------------------------//
531  // Z_eq_b //
532  //----------------------------//
533 
534 double Star_rot::z_eqb() const {
535 
536  if (p_z_eqb == 0x0) { // a new computation is required
537 
538  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
539  const Mg3d* mg = mp.get_mg() ;
540  assert(mg->get_type_t() == SYM) ;
541  int l_b = nzet - 1 ;
542  int i_b = mg->get_nr(l_b) - 1 ;
543  int j_b = mg->get_nt(l_b) - 1 ;
544  int k_b = 0 ;
545 
546  double u_eq = uuu.val_grid_point(l_b, k_b, j_b, i_b) ;
547  double n_eq = nn.val_grid_point(l_b, k_b, j_b, i_b) ;
548  double nphi_eq = nphi.val_grid_point(l_b, k_b, j_b, i_b) ;
549 
550  p_z_eqb = new double( sqrt((1.+u_eq)/(1.-u_eq))
551  / (n_eq - nphi_eq * r_circ() )
552  - 1. ) ;
553 
554  }
555 
556  return *p_z_eqb ;
557 
558 }
559 
560 
561  //----------------------------//
562  // Z_pole //
563  //----------------------------//
564 
565 double Star_rot::z_pole() const {
566 
567  if (p_z_pole == 0x0) { // a new computation is required
568 
569  double n_pole = nn.val_point(ray_pole(), 0., 0.) ;
570 
571  p_z_pole = new double( 1. / n_pole - 1. ) ;
572 
573  }
574 
575  return *p_z_pole ;
576 
577 }
578 
579 
580  //----------------------------//
581  // Quadrupole moment //
582  // Uses Eq.(11) of Pappas & //
583  // Apostolatos, PRL 108 2012 //
584  //----------------------------//
585 
586 double Star_rot::mom_quad() const {
587 
588  using namespace Unites ;
589  if (p_mom_quad == 0x0) { // a new computation is required
590 
591  // Source for of the Poisson equation for nu
592  // -----------------------------------------
593 
594  Scalar source(mp) ;
595 
596  if (relativistic) {
597  Scalar bet = log(bbb) ;
598  bet.std_spectral_base() ;
599 
600  Vector dlogn = logn.derive_cov( mp.flat_met_spher() ) ;
601  Vector dlogn_bet = dlogn + bet.derive_cov( mp.flat_met_spher() ) ;
602 
603  source = qpig * a_car *( ener_euler + s_euler ) + ak_car
604  - dlogn(1)*dlogn_bet(1) - dlogn(2)*dlogn_bet(2) - dlogn(3)*dlogn_bet(3) ;
605  }
606  else {
607  source = qpig * nbar ;
608  }
609 
610  source.std_spectral_base() ;
611 
612  // Multiplication by -r^2 P_2(cos(theta))
613  // [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
614  // ------------------------------------------------------------------
615 
616  // Multiplication by r^2 :
617  // ----------------------
618  source.mult_r() ;
619  source.mult_r() ;
620  if (source.check_dzpuis(2)) {
621  source.inc_dzpuis(2) ;
622  }
623 
624  // Muliplication by cos^2(theta) :
625  // -----------------------------
626  Scalar temp = source ;
627 
628  // What follows is valid only for a mapping of class Map_radial :
629  assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ;
630 
631  if (temp.get_etat() == ETATQCQ) {
632  Valeur& vtemp = temp.set_spectral_va() ;
633  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
634  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
635  }
636 
637  // Muliplication by -P_2(cos(theta)) :
638  // ----------------------------------
639  source = 0.5 * source - 1.5 * temp ;
640 
641  // First part
642  // ------------
643 
644  p_mom_quad = new double( -source.integrale() / qpig ) ;
645 
646  // Correction terms from (3.37) of Friedman & Stergioulas' book
647  if (relativistic) {
648 
649  Scalar correct(mp) ;
650  correct = a_car * bbb * nn * press ;
651  correct.mult_rsint() ;
652  correct.std_spectral_base() ;
653 
654  double Bo = -32. * correct.integrale() / qpig ;
655 
656  *p_mom_quad -= 4./3. * (1./4. + Bo/(mass_g() * mass_g()))
657  * pow(mass_g(), 3) * ggrav * ggrav ;
658  }
659 
660  }
661 
662  return *p_mom_quad ;
663 
664 }
665 
666 
667  //----------------------------//
668  // Surface gravity //
669  //----------------------------//
670 
671 const Tbl& Star_rot::surf_grav() const {
672 
673  using namespace Unites ;
674  if (p_surf_grav == 0x0) { // a new computation is required
675 
676  const Mg3d* mg = mp.get_mg() ;
677  assert(mg->get_type_t() == SYM) ;
678 
679  // Index of the point at phi=0 at the surface of the star:
680  int l_b = nzet - 1 ;
681  int i_b = mg->get_nr(l_b) - 1 ;
682  int j_b = mg->get_nt(l_b) ; // number of theta points
683 
684  Scalar loggam = log(gam_euler) ;
685  loggam.std_spectral_base() ;
686 
688  - loggam.derive_cov( mp.flat_met_spher() );
689 
690  Scalar g = contract( a, 0, a.up_down(gamma), 0 ) ;
691 
692  p_surf_grav = new Tbl(j_b) ;
694 
695  int j;
696  for(j = 0; j < j_b; j++)
697  p_surf_grav->set(j) = sqrt(g.val_grid_point(l_b, 0, j, i_b)) ;
698 
699  }
700 
701  return *p_surf_grav ;
702 
703 }
704 
705 }
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:333
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:146
virtual double mean_radius() const
Mean star radius from the area .
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:116
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:479
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
double * p_r_circ
Circumferential radius (equator)
Definition: star_rot.h:248
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:106
Scalar tnphi
Component of the shift vector.
Definition: star_rot.h:130
double * p_mom_quad
Quadrupole moment.
Definition: star_rot.h:255
Scalar bbb
Metric factor B.
Definition: star_rot.h:119
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:783
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:393
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:64
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:502
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:113
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:560
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
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:643
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:245
Scalar nphi
Metric coefficient .
Definition: star_rot.h:125
double * p_z_eqf
Forward redshift factor at equator.
Definition: star_rot.h:252
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:246
double * p_area
Integrated surface area.
Definition: star_rot.h:251
double * p_r_circ_merid
Circumferential radius (meridian)
Definition: star_rot.h:249
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:604
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:563
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:254
Base class for pure radial mappings.
Definition: map.h:1557
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:330
double * p_grv3
Error on the virial identity GRV3.
Definition: star_rot.h:247
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:1570
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:244
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:122
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:469
double * p_aplat
Flatening r_pole/r_eq.
Definition: star_rot.h:250
Multi-domain grid.
Definition: grilles.h:279
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:145
Scalar nn
Lapse function N .
Definition: star.h:225
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:133
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:113
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:474
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:610
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:198
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:390
Tbl * p_surf_grav
Surface gravity (along the theta direction)
Definition: star_rot.h:263
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:324
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:607
double * p_z_eqb
Backward redshift factor at equator.
Definition: star_rot.h:253