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