LORENE
et_bfrot_global.C
1 /*
2  * Copyright (c) 2001 Jerome Novak
3  *
4  * This file is part of LORENE.
5  *
6  * LORENE is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * LORENE is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with LORENE; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  *
20  */
21 
22 
23 
24 
25 /*
26  * $Id: et_bfrot_global.C,v 1.19 2017/10/06 12:36:34 a_sourie Exp $
27  * $Log: et_bfrot_global.C,v $
28  * Revision 1.19 2017/10/06 12:36:34 a_sourie
29  * Cleaning of tabulated 2-fluid EoS class + superfluid rotating star model.
30  *
31  * Revision 1.18 2016/12/05 16:17:52 j_novak
32  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
33  *
34  * Revision 1.17 2015/06/10 14:39:17 a_sourie
35  * New class Eos_bf_tabul for tabulated 2-fluid EoSs and associated functions for the computation of rotating stars with such EoSs.
36  *
37  * Revision 1.16 2014/10/13 08:52:54 j_novak
38  * Lorene classes and functions now belong to the namespace Lorene.
39  *
40  * Revision 1.15 2014/10/06 15:13:07 j_novak
41  * Modified #include directives to use c++ syntax.
42  *
43  * Revision 1.14 2004/09/03 13:55:07 j_novak
44  * Use of enerps_euler instead of ener_euler in the calculation of mom_angu()
45  *
46  * Revision 1.13 2004/03/25 10:29:03 j_novak
47  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
48  *
49  * Revision 1.12 2003/11/20 14:01:26 r_prix
50  * changed member names to better conform to Lorene coding standards:
51  * J_euler -> j_euler, EpS_euler -> enerps_euler, Delta_car -> delta_car
52  *
53  * Revision 1.11 2003/11/18 18:38:11 r_prix
54  * use of new member EpS_euler: matter sources in equilibrium() and global quantities
55  * no longer distinguish Newtonian/relativistic, as all terms should have the right limit...
56  *
57  * Revision 1.10 2003/11/13 12:13:15 r_prix
58  * *) removed all uses of Etoile-type specific quantities: u_euler, press
59  * and exclusively use ener_euler, s_euler, sphph_euler, J_euler instead
60  * *) this also fixes a bug in previous expression for mass_g() in 2-fluid case
61  *
62  * NOTE: some current Newtonian expressions are still completely broken..
63  *
64  * Revision 1.9 2003/09/17 08:27:50 j_novak
65  * New methods: mass_b1() and mass_b2().
66  *
67  * Revision 1.8 2003/02/07 10:47:43 j_novak
68  * The possibility of having gamma5 xor gamma6 =0 has been introduced for
69  * tests. The corresponding parameter files have been added.
70  *
71  * Revision 1.7 2002/10/16 14:36:36 j_novak
72  * Reorganization of #include instructions of standard C++, in order to
73  * use experimental version 3 of gcc.
74  *
75  * Revision 1.6 2002/10/14 14:20:08 j_novak
76  * Error corrected for angu_mom()
77  *
78  * Revision 1.5 2002/05/10 09:26:52 j_novak
79  * Added new class Et_rot_mag for magnetized rotating neutron stars (under development)
80  *
81  * Revision 1.4 2002/04/05 09:09:36 j_novak
82  * The inversion of the EOS for 2-fluids polytrope has been modified.
83  * Some errors in the determination of the surface were corrected.
84  *
85  * Revision 1.3 2002/01/08 14:43:53 j_novak
86  * better determination of surfaces for 2-fluid stars
87  *
88  * Revision 1.2 2002/01/03 15:30:28 j_novak
89  * Some comments modified.
90  *
91  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
92  * LORENE
93  *
94  * Revision 1.4 2001/08/31 15:07:12 novak
95  * Retour a la version 1.2, sans la routine prolonge_c1
96  *
97  * Revision 1.2 2001/08/28 15:32:17 novak
98  * The surface is now defined by the baryonic density
99  *
100  * Revision 1.1 2001/06/22 15:39:46 novak
101  * Initial revision
102  *
103  *
104  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_global.C,v 1.19 2017/10/06 12:36:34 a_sourie Exp $
105  *
106  */
107 
108 
109 // Headers C
110 #include <cstdlib>
111 #include <cmath>
112 
113 // Headers Lorene
114 #include "et_rot_bifluid.h"
115 #include "unites.h"
116 
117 //--------------------------//
118 // Baryon mass //
119 //--------------------------//
120 
121 namespace Lorene {
122 double Et_rot_bifluid::mass_b1() const {
123 
124  if (p_mass_b1 == 0x0) { // a new computation is required
125 
126  assert(nbar.get_etat() == ETATQCQ);
127  Cmp dens1 = a_car() * bbb() * gam_euler() * eos.get_m1()* nbar();
128  dens1.std_base_scal() ;
129 
130  p_mass_b1 = new double( dens1.integrale() ) ;
131 
132  }
133 
134  return *p_mass_b1 ;
135 
136 }
137 
138 double Et_rot_bifluid::mass_b2() const {
139 
140  if (p_mass_b2 == 0x0) { // a new computation is required
141  assert(nbar2.get_etat() == ETATQCQ);
142 
143  Cmp dens2 = a_car() * bbb() * gam_euler2() * eos.get_m2() * nbar2();
144  dens2.std_base_scal() ;
145 
146  p_mass_b2 = new double( dens2.integrale() ) ;
147 
148  }
149 
150  return *p_mass_b2 ;
151 
152 }
153 
154 double Et_rot_bifluid::mass_b() const {
155 
156  if (p_mass_b == 0x0)
157  p_mass_b = new double(mass_b1() + mass_b2() ) ;
158 
159  return *p_mass_b;
160 
161 }
162 
163 
164 //----------------------------//
165 // Gravitational mass //
166 //----------------------------//
167 
168 double Et_rot_bifluid::mass_g() const {
169 
170  if (p_mass_g == 0x0) { // a new computation is required
171 
172  Tenseur vtmp (j_euler);
173  vtmp.change_triad( mp.get_bvect_spher() ); // change to spherical triad
174 
175  Tenseur tJphi (mp);
176  tJphi = vtmp(2); // get the phi tetrad-component: J^ph r sin(th)
177 
178  // relativistic AND Newtonian: enerps_euler has right limit and tnphi->0
179  Tenseur source = nnn * enerps_euler + 2 * b_car * tnphi * tJphi;
180  source = a_car * bbb * source ;
181 
182  source.set_std_base() ;
183 
184  p_mass_g = new double( source().integrale() ) ;
185 
186  }
187 
188  return *p_mass_g ;
189 
190 }
191 
192 //----------------------------//
193 // Angular momentum //
194 //----------------------------//
195 
196 double Et_rot_bifluid::angu_mom() const {
197 
198  if (p_angu_mom == 0x0) { // a new computation is required
199 
200  Cmp dens(mp) ;
201 
202  // this should work for both relativistic AND Newtonian,
203  // provided j_euler has the right limit...
204  Tenseur tmp = j_euler;
205  tmp.change_triad( mp.get_bvect_spher() ) ;
206  dens = tmp(2) ;
207  dens.mult_rsint() ;
208 
209  dens = a_car() * b_car()* bbb() * dens ;
210 
211  dens.std_base_scal() ;
212 
213  p_angu_mom = new double( dens.integrale() ) ;
214 
215  }
216 
217  return *p_angu_mom ;
218 
219 }
220 
221 
223 
224  if (p_angu_mom_1 == 0x0) { // a new computation is required
225 
226  Cmp dens(mp) ;
227 
228  // this should work for both relativistic AND Newtonian,
229  // provided j_euler has the right limit...
230  Tenseur tmp = j_euler1;
231  tmp.change_triad( mp.get_bvect_spher() ) ;
232  dens = tmp(2) ;
233  dens.mult_rsint() ;
234 
235  dens = a_car() * b_car()* bbb() * dens ;
236 
237  dens.std_base_scal() ;
238 
239  p_angu_mom_1 = new double( dens.integrale() ) ;
240 
241  }
242 
243  return *p_angu_mom_1 ;
244 
245 }
246 
248 
249  if (p_angu_mom_2 == 0x0) { // a new computation is required
250 
251  Cmp dens(mp) ;
252 
253  // this should work for both relativistic AND Newtonian,
254  // provided j_euler has the right limit...
255  Tenseur tmp = j_euler2;
256  tmp.change_triad( mp.get_bvect_spher() ) ;
257  dens = tmp(2) ;
258  dens.mult_rsint() ;
259 
260  dens = a_car() * b_car()* bbb() * dens ;
261 
262  dens.std_base_scal() ;
263 
264  p_angu_mom_2 = new double( dens.integrale() ) ;
265 
266  }
267 
268  return *p_angu_mom_2 ;
269 
270 }
271 
272 //--------------------------//
273 // Fluid couplings //
274 //--------------------------//
275 
276 /* The following quantities are defined in Eqs. (C5) - (C7) of
277 Sourie et al., MNRAS 464 (2017).
278 Stricly speaking, they only make sense in the slow-rotation
279 approximation and to first order in the lag
280 */
281 
283 
284  if (p_coupling_mominert_1 == 0x0) { // a new computation is required
285 
286  Cmp dens(mp) ;
287 
288  Tenseur mu_1 = eos.get_m1() * exp( unsurc2* ent);
289  mu_1.set_std_base() ;
290 
291  //dens = gam_euler() * gam_euler() * nbar() * mu_1() * a_car() * b_car() * bbb() / nnn() ;
292  dens = nbar() * mu_1() * a_car() * b_car() * bbb() / nnn() ;
293  dens.std_base_scal() ;
294  dens.mult_rsint() ;
295  dens.std_base_scal() ;
296  dens.mult_rsint() ;
297  dens.std_base_scal() ;
298 
299  p_coupling_mominert_1 = new double( dens.integrale() );
300 
301  }
302 
303  return *p_coupling_mominert_1 ;
304 
305 }
306 
307 
308 double Et_rot_bifluid::coupling_mominert_2() const {
309 
310  if (p_coupling_mominert_2 == 0x0) { // a new computation is required
311 
312  Cmp dens(mp) ;
313 
314  Tenseur mu_2 = eos.get_m2() * exp(unsurc2 * ent2);
315  mu_2.set_std_base() ;
316 
317  //dens = gam_euler2() * gam_euler2() * nbar2() * mu_2() * a_car() * b_car() * bbb() / nnn() ;
318  dens = nbar2() * mu_2() * a_car() * b_car() * bbb() / nnn() ;
319  dens.std_base_scal() ;
320  dens.mult_rsint() ;
321  dens.std_base_scal() ;
322  dens.mult_rsint() ;
323  dens.std_base_scal() ;
324 
325  p_coupling_mominert_2 = new double( dens.integrale() );
326 
327  }
328 
329  return *p_coupling_mominert_2 ;
330 
331 }
332 
333 
334 double Et_rot_bifluid::coupling_entr() const {
335 
336  if (p_coupling_entr == 0x0) { // a new computation is required
337 
338  Cmp dens(mp) ;
339 
340  Tenseur gam_delta = 1./sqrt(1.-unsurc2 * delta_car);
341 
342  dens = 2. * alpha_eos() * a_car() * b_car() * bbb() / nnn() ;
343  dens.std_base_scal() ;
344  dens.mult_rsint() ;
345  dens.std_base_scal() ;
346  dens.mult_rsint() ;
347  dens.std_base_scal() ;
348 
349  p_coupling_entr = new double( dens.integrale() );
350 
351  }
352 
353  return *p_coupling_entr ;
354 
355 }
356 
357 double Et_rot_bifluid::coupling_LT_1() const {
358 
359  if (p_coupling_LT_1 == 0x0) { // a new computation is required
360 
361  Cmp dens(mp) ;
362 
363  Tenseur mu_1 = eos.get_m1() * exp( unsurc2* ent); // tester avec Knn et Knp
364  mu_1.set_std_base() ;
365 
366  // dens = gam_euler() * gam_euler() * nbar() * mu_1() * a_car() * b_car() * bbb() * nphi()/ nnn() ;
367  dens = nbar() * mu_1() * a_car() * b_car() * bbb() * nphi()/ nnn() ;
368  dens.std_base_scal() ;
369  dens.mult_rsint() ;
370  dens.std_base_scal() ;
371  dens.mult_rsint() ;
372  dens.std_base_scal() ;
373 
374  p_coupling_LT_1 = new double( dens.integrale() );
375 
376  }
377 
378  return *p_coupling_LT_1 ;
379 
380 }
381 
382 
383 double Et_rot_bifluid::coupling_LT_2() const {
384 
385  if (p_coupling_LT_2 == 0x0) { // a new computation is required
386 
387  Cmp dens(mp) ;
388 
389  Tenseur mu_2 = eos.get_m2() * exp(unsurc2 * ent2); // tester avec Knn et Knp
390  mu_2.set_std_base() ;
391 
392  // dens = gam_euler2() * gam_euler2() * nbar2() * mu_2() * a_car() * b_car() * bbb() * nphi()/ nnn() ;
393  dens = nbar2() * mu_2() * a_car() * b_car() * bbb() * nphi()/ nnn() ;
394  dens.std_base_scal() ;
395  dens.mult_rsint() ;
396  dens.std_base_scal() ;
397  dens.mult_rsint() ;
398  dens.std_base_scal() ;
399 
400  p_coupling_LT_2 = new double( dens.integrale() );
401 
402  }
403 
404  return *p_coupling_LT_2 ;
405 
406 }
407 
408 //----------------------------//
409 // GRV2 //
410 //----------------------------//
411 
412 double Et_rot_bifluid::grv2() const {
413 
414  using namespace Unites ;
415 
416  if (p_grv2 == 0x0) { // a new computation is required
417 
418  Tenseur sou_m(mp);
419  // should work for both relativistic AND Newtonian, provided sphph_euler has right limit..
420  sou_m = 2 * qpig * a_car * sphph_euler ;
421 
423 
424  p_grv2 = new double( double(1) - lambda_grv2(sou_m(), sou_q()) ) ;
425 
426  }
427 
428  return *p_grv2 ;
429 
430 }
431 
432 
433 //----------------------------//
434 // GRV3 //
435 //----------------------------//
436 
437 double Et_rot_bifluid::grv3(ostream* ost) const {
438 
439  using namespace Unites ;
440 
441  if (p_grv3 == 0x0) { // a new computation is required
442 
443  Tenseur source(mp) ;
444 
445  // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
446  // ------------------ Class. Quantum Grav. 11, 443 (1994)]
447 
448  if (relativistic) {
449  Tenseur alpha = dzeta - logn ;
450  Tenseur beta = log( bbb ) ;
451  beta.set_std_base() ;
452 
454  + 0.5 * flat_scalar_prod(alpha.gradient_spher(), beta.gradient_spher() ) ;
455 
456  Cmp aa = alpha() - 0.5 * beta() ;
457  Cmp daadt = aa.srdsdt() ; // 1/r d/dth
458 
459  // What follows is valid only for a mapping of class Map_radial :
460  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
461  if (mpr == 0x0) {
462  cout << "Et_rot_bifluid::grv3: the mapping does not belong"
463  << " to the class Map_radial !" << endl ;
464  abort() ;
465  }
466 
467  // Computation of 1/tan(theta) * 1/r daa/dtheta
468  if (daadt.get_etat() == ETATQCQ) {
469  Valeur& vdaadt = daadt.va ;
470  vdaadt = vdaadt.ssint() ; // division by sin(theta)
471  vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
472  }
473 
474  Cmp temp = aa.dsdr() + daadt ;
475  temp = ( bbb() - a_car()/bbb() ) * temp ;
476  temp.std_base_scal() ;
477 
478  // Division by r
479  Valeur& vtemp = temp.va ;
480  vtemp = vtemp.sx() ; // division by xi in the nucleus
481  // Id in the shells
482  // division by xi-1 in the ZEC
483  vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
484  // by 1/r in the shells
485  // by r(xi-1) in the ZEC
486 
487  // In the ZEC, a multiplication by r has been performed instead
488  // of the division:
489  temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;
490 
491  source = bbb() * source() + 0.5 * temp ;
492 
493  }
494  else{
495  source = - 0.5 * flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() ) ;
496  }
497 
498  source.set_std_base() ;
499 
500  double int_grav = source().integrale() ;
501 
502  // Matter term
503  // -----------
504 
505  // should work for relativistic AND Newtonian, provided s_euler has the right limit...
506  source = qpig * a_car * bbb * s_euler ;
507 
508  source.set_std_base() ;
509 
510  double int_mat = source().integrale() ;
511 
512  // Virial error
513  // ------------
514  if (ost != 0x0) {
515  *ost << "Et_rot_bifluid::grv3 : gravitational term : " << int_grav
516  << endl ;
517  *ost << "Et_rot_bifluid::grv3 : matter term : " << int_mat
518  << endl ;
519  }
520 
521  p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;
522 
523  }
524 
525  return *p_grv3 ;
526 
527 }
528 
529 
530 //----------------------------//
531 // R_circ //
532 //----------------------------//
533 
534 double Et_rot_bifluid::r_circ2() const {
535 
536  if (p_r_circ2 == 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  p_r_circ2 = new double( bbb()(l_b, k_b, j_b, i_b) * ray_eq2() ) ;
547 
548  }
549 
550  return *p_r_circ2 ;
551 
552 }
553 
554 
555 //----------------------------//
556 // Flattening //
557 //----------------------------//
558 
559 double Et_rot_bifluid::aplat2() const {
560 
561  if (p_aplat2 == 0x0) { // a new computation is required
562 
563  p_aplat2 = new double( ray_pole2() / ray_eq2() ) ;
564 
565  }
566 
567  return *p_aplat2 ;
568 
569 }
570 
571 
572 
573 //----------------------------//
574 // Surface area //
575 //----------------------------//
576 
577  double Et_rot_bifluid::area2() const {
578 
579  if (p_area2 == 0x0) { // a new computation is required
580  const Mg3d& mg = *(mp.get_mg()) ;
581  int np = mg.get_np(0) ;
582  int nt = mg.get_nt(0) ;
583  assert(np == 1) ; //Only valid for axisymmetric configurations
584 
585  const Map_radial* mp_rad = dynamic_cast<const Map_radial*>(&mp) ;
586  assert(mp_rad != 0x0) ;
587 
588  Valeur va_r(mg.get_angu()) ;
589  va_r.annule_hard() ;
590  Itbl lsurf2 = l_surf2() ;
591  Tbl xisurf2 = xi_surf2() ;
592 
593  for (int k=0; k<np; k++) {
594  for (int j=0; j<nt; j++) {
595  int l_star2 = lsurf2(k, j) ;
596  double xi_star2 = xisurf2(k, j) ;
597 
598  va_r.set(0, k, j, 0) = mp_rad->val_r_jk(l_star2, xi_star2, j, k) ;
599  }
600  }
601  va_r.std_base_scal() ;
602 
603  Valeur integ(mg.get_angu()) ;
604  Valeur dr = va_r.dsdt() ;
605  integ = sqrt(va_r*va_r + dr*dr) ;
606  Cmp aaaa = get_a_car()() ;
607  Valeur a2 = aaaa.va ; a2.std_base_scal() ;
608  Cmp bbbb = get_bbb()() ;
609  Valeur b = bbbb.va ; b.std_base_scal() ;
610  for (int k=0; k<np; k++) {
611  for (int j=0; j<nt; j++) {
612  integ.set(0, k, j, 0) *= sqrt(a2.val_point_jk(lsurf2(k, j), xisurf2(k, j), j, k))
613  * b.val_point_jk(lsurf2(k, j), xisurf2(k, j), j, k) * va_r(0, k, j, 0) ;
614  }
615  }
616  integ.std_base_scal() ;
617  Valeur integ2 = integ.mult_st() ;
618  double surftot = 0. ;
619  for (int j=0; j<nt; j++) {
620  surftot += (*integ2.c_cf)(0, 0, j, 0) / double(2*j+1) ;
621  }
622 
623  p_area2 = new double( 4*M_PI*surftot) ;
624 
625  }
626 
627  return *p_area2 ;
628 
629 }
630 
632 
633  return sqrt(area2()/(4*M_PI)) ;
634 
635  }
636 
637 
638 
639 
640 //----------------------------//
641 // Quadrupole moment //
642 //----------------------------//
643 
644 double Et_rot_bifluid::mom_quad() const {
645 
646  using namespace Unites ;
647 
648  if (p_mom_quad == 0x0) { // a new computation is required
649 
650  double b = mom_quad_Bo() / ( mass_g() * mass_g() ) ;
651  p_mom_quad = new double( mom_quad_old() - 4./3. * ( 1./4. + b ) * pow(mass_g(), 3) * ggrav * ggrav ) ;
652 
653  }
654 
655  return *p_mom_quad ;
656 
657 }
658 
659 
661 
662  using namespace Unites ;
663 
664  if (p_mom_quad_Bo == 0x0) { // a new computation is required
665 
666  Cmp dens(mp) ;
667 
668  dens = press() ;
669  dens = a_car() * bbb() * nnn() * dens ;
670  dens.mult_rsint() ;
671  dens.std_base_scal() ;
672 
673  p_mom_quad_Bo = new double( - 32. * dens.integrale() / qpig ) ;
674 
675  }
676 
677  return *p_mom_quad_Bo ;
678 
679 }
680 
681 
682 
684 
685  using namespace Unites ;
686 
687  if (p_mom_quad_old == 0x0) { // a new computation is required
688 
689  // Source for of the Poisson equation for nu
690  // -----------------------------------------
691 
692  Tenseur source(mp) ;
693 
694  if (relativistic) {
695  Tenseur beta = log(bbb) ;
696  beta.set_std_base() ;
697  source = qpig * a_car * enerps_euler
699  logn.gradient_spher() + beta.gradient_spher()) ;
700  }
701  else {
702  source = qpig * (nbar + nbar2);
703  }
704  source.set_std_base() ;
705 
706  // Multiplication by -r^2 P_2(cos(theta))
707  // [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
708  // ------------------------------------------------------------------
709 
710  // Multiplication by r^2 :
711  // ----------------------
712  Cmp& csource = source.set() ;
713  csource.mult_r() ;
714  csource.mult_r() ;
715  if (csource.check_dzpuis(2)) {
716  csource.inc2_dzpuis() ;
717  }
718 
719  // Muliplication by cos^2(theta) :
720  // -----------------------------
721  Cmp temp = csource ;
722 
723  // What follows is valid only for a mapping of class Map_radial :
724  assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ;
725 
726  if (temp.get_etat() == ETATQCQ) {
727  Valeur& vtemp = temp.va ;
728  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
729  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
730  }
731 
732  // Muliplication by -P_2(cos(theta)) :
733  // ----------------------------------
734  source = 0.5 * source() - 1.5 * temp ;
735 
736  // Final result
737  // ------------
738 
739  p_mom_quad_old = new double(- source().integrale() / qpig ) ;
740 
741  }
742 
743  return *p_mom_quad_old ;
744 
745 }
746 
747 
748 //--------------------------//
749 // Stellar surface //
750 //--------------------------//
751 
752 const Itbl& Et_rot_bifluid::l_surf() const {
753 
754  if (p_l_surf == 0x0) { // a new computation is required
755 
756  assert(p_xi_surf == 0x0) ; // consistency check
757 
758  int np = mp.get_mg()->get_np(0) ;
759  int nt = mp.get_mg()->get_nt(0) ;
760 
761  p_l_surf = new Itbl(np, nt) ;
762  p_xi_surf = new Tbl(np, nt) ;
763 
764  double nb0 = 0 ; // definition of the surface
765  double precis = 1.e-15 ;
766  int nitermax = 10000000 ;
767  int niter ;
768 
769  // Cmp defining the surface of the star (via the density fields)
770  //
771  Cmp surf(mp) ;
772  surf = -0.2*nbar()(0,0,0,0) ;
773  surf.annule(0, nzet-1) ;
774  surf += nbar() ; ;
775  surf = prolonge_c1(surf, nzet) ;
776 
777  (surf.va).equipot(nb0, nzet, precis, nitermax, niter, *p_l_surf,
778  *p_xi_surf) ;
779 
780  }
781 
782  return *p_l_surf ;
783 
784 }
786 
787  if (p_l_surf2 == 0x0) { // a new computation is required
788 
789  assert(p_xi_surf2 == 0x0) ; // consistency check
790 
791  int np = mp.get_mg()->get_np(0) ;
792  int nt = mp.get_mg()->get_nt(0) ;
793 
794  p_l_surf2 = new Itbl(np, nt) ;
795  p_xi_surf2 = new Tbl(np, nt) ;
796 
797  double nb0 = 0 ; // definition of the surface
798  double precis = 1.e-15 ;
799  int nitermax = 10000000 ;
800  int niter ;
801 
802  // Cmp defining the surface of the star (via the density fields)
803  //
804  Cmp surf2(mp) ;
805  surf2 = -0.2*nbar2()(0,0,0,0) ;
806  surf2.annule(0, nzet-1) ;
807  surf2 += nbar2() ; ;
808  surf2 = prolonge_c1(surf2, nzet) ;
809 
810  (surf2.va).equipot(nb0, nzet, precis, nitermax, niter, *p_l_surf2,
811  *p_xi_surf2) ;
812 
813  }
814 
815  return *p_l_surf2 ;
816 
817 }
818 
820 
821  if (p_xi_surf2 == 0x0) { // a new computation is required
822 
823  assert(p_l_surf2 == 0x0) ; // consistency check
824 
825  l_surf2() ; // the computation is done by l_surf2()
826 
827  }
828 
829  return *p_xi_surf2 ;
830 
831 }
832 
833 //--------------------------//
834 // Coordinate radii //
835 //--------------------------//
836 
837 double Et_rot_bifluid::ray_eq2() const {
838 
839  if (p_ray_eq2 == 0x0) { // a new computation is required
840 
841  const Mg3d& mg = *(mp.get_mg()) ;
842 
843  int type_t = mg.get_type_t() ;
844  int nt = mg.get_nt(0) ;
845 
846  if ( type_t == SYM ) {
847  assert( ( mg.get_type_p() == SYM) || (mg.get_type_p() == NONSYM) ) ;
848  int k = 0 ;
849  int j = nt-1 ;
850  int l = l_surf2()(k, j) ;
851  double xi = xi_surf2()(k, j) ;
852  double theta = M_PI / 2 ;
853  double phi = 0 ;
854 
855  p_ray_eq2 = new double( mp.val_r(l, xi, theta, phi) ) ;
856 
857  }
858  else {
859  cout << "Et_rot_bifluid::ray_eq2 : the case type_t = " << type_t
860  << " is not contemplated yet !" << endl ;
861  abort() ;
862  }
863 
864  }
865 
866  return *p_ray_eq2 ;
867 
868 }
869 
870 
872 
873  if (p_ray_eq2_pis2 == 0x0) { // a new computation is required
874 
875  const Mg3d& mg = *(mp.get_mg()) ;
876 
877  int type_t = mg.get_type_t() ;
878  int type_p = mg.get_type_p() ;
879  int nt = mg.get_nt(0) ;
880  int np = mg.get_np(0) ;
881 
882  if ( type_t == SYM ) {
883 
884  int j = nt-1 ;
885  double theta = M_PI / 2 ;
886  double phi = M_PI / 2 ;
887 
888  switch (type_p) {
889 
890  case SYM : {
891  int k = np / 2 ;
892  int l = l_surf2()(k, j) ;
893  double xi = xi_surf2()(k, j) ;
894  p_ray_eq2_pis2 = new double( mp.val_r(l, xi, theta, phi) ) ;
895  break ;
896  }
897 
898  case NONSYM : {
899  assert( np % 4 == 0 ) ;
900  int k = np / 4 ;
901  int l = l_surf2()(k, j) ;
902  double xi = xi_surf2()(k, j) ;
903  p_ray_eq2_pis2 = new double( mp.val_r(l, xi, theta, phi) ) ;
904  break ;
905  }
906 
907  default : {
908  cout << "Et_rot_bifluid::ray_eq2_pis2 : the case type_p = "
909  << type_p << " is not contemplated yet !" << endl ;
910  abort() ;
911  }
912  }
913 
914  }
915  else {
916  cout << "Et_rot_bifluid::ray_eq2_pis2 : the case type_t = " << type_t
917  << " is not contemplated yet !" << endl ;
918  abort() ;
919  }
920 
921  }
922 
923  return *p_ray_eq2_pis2 ;
924 
925 }
926 
927 
929 
930  if (p_ray_eq2_pi == 0x0) { // a new computation is required
931 
932  const Mg3d& mg = *(mp.get_mg()) ;
933 
934  int type_t = mg.get_type_t() ;
935  int type_p = mg.get_type_p() ;
936  int nt = mg.get_nt(0) ;
937  int np = mg.get_np(0) ;
938 
939  if ( type_t == SYM ) {
940 
941  switch (type_p) {
942 
943  case SYM : {
944  p_ray_eq2_pi = new double( ray_eq2() ) ;
945  break ;
946  }
947 
948  case NONSYM : {
949  int k = np / 2 ;
950  int j = nt-1 ;
951  int l = l_surf2()(k, j) ;
952  double xi = xi_surf2()(k, j) ;
953  double theta = M_PI / 2 ;
954  double phi = M_PI ;
955 
956  p_ray_eq2_pi = new double( mp.val_r(l, xi, theta, phi) ) ;
957  break ;
958  }
959 
960  default : {
961 
962  cout << "Et_rot_bifluid::ray_eq2_pi : the case type_t = " << type_t
963  << " and type_p = " << type_p << endl ;
964  cout << " is not contemplated yet !" << endl ;
965  abort() ;
966  break ;
967  }
968  }
969  }
970 
971  }
972 
973  return *p_ray_eq2_pi ;
974 
975 }
976 
978 
979  if (p_ray_pole2 == 0x0) { // a new computation is required
980 
981  assert( ((mp.get_mg())->get_type_t() == SYM)
982  || ((mp.get_mg())->get_type_t() == NONSYM) ) ;
983 
984  int k = 0 ;
985  int j = 0 ;
986  int l = l_surf2()(k, j) ;
987  double xi = xi_surf2()(k, j) ;
988  double theta = 0 ;
989  double phi = 0 ;
990 
991  p_ray_pole2 = new double( mp.val_r(l, xi, theta, phi) ) ;
992 
993  }
994 
995  return *p_ray_pole2 ;
996 
997 }
998 
999 
1000 
1001 }
int get_type_p() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the transformatio...
Definition: grilles.h:512
const Cmp & dsdr() const
Returns of *this .
Definition: cmp_deriv.C:87
double get_m2() const
Return the individual particule mass .
Definition: eos_bifluid.h:269
double * p_mass_b2
Baryon mass of fluid 2.
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_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
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
virtual double mass_g() const
Gravitational mass.
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
Tenseur j_euler2
To compute .
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition: tenseur.C:1564
Tenseur j_euler1
To compute .
virtual double mass_b() const
Total Baryon mass.
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 mass_b1() const
Baryon mass of fluid 1.
void annule(int l)
Sets the Cmp to zero in a given domain.
Definition: cmp.C:351
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:801
double * p_r_circ2
Circumferential radius of fluid no.2.
Tbl * p_xi_surf2
Description of the surface of fluid 2: 2-D Tbl containing the values of the radial coordinate on the...
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
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 * p_angu_mom_1
Angular momentum of fluid 1.
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
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
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition: etoile.h:445
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
double * p_ray_eq2
Coordinate radius at , .
Tenseur press
Fluid pressure.
Definition: etoile.h:464
virtual double mom_quad_Bo() const
Part of the quadrupole moment.
const Itbl & l_surf2() const
Description of the surface of fluid 2: returns a 2-D Itbl containing the values of the domain index ...
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition: valeur.C:827
Tenseur j_euler
Total angular momentum (flat-space!) 3-vector , which is related to of the "3+1" decomposition...
void mult_r()
Multiplication by r everywhere.
Definition: cmp_r_manip.C:94
double * p_mass_b1
Baryon mass of fluid 1.
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
Tenseur gam_euler2
Lorentz factor between the fluid 2 and Eulerian observers.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
double * p_grv3
Error on the virial identity GRV3.
Definition: etoile.h:1637
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:684
double ray_eq2_pi() const
Coordinate radius for fluid 2 at , [r_unit].
const Eos_bifluid & eos
Equation of state for two-fluids model.
const Valeur & ssint() const
Returns of *this.
Definition: valeur_ssint.C:115
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
double * p_coupling_mominert_1
Quantities used to describe the different couplings between the fluids.
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:474
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Base class for pure radial mappings.
Definition: map.h:1557
virtual double angu_mom_2() const
Angular momentum of fluid 2.
void mult_rsint()
Multiplication by .
Definition: cmp_r_manip.C:119
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
double * p_area2
Surface area of fluid no.2.
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
const Tbl & xi_surf2() const
Description of the surface of fluid 2: returns a 2-D Tbl containing the values of the radial coordin...
double ray_eq2_pis2() const
Coordinate radius for fluid 2 at , [r_unit].
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.
Tenseur enerps_euler
the combination : useful because in the Newtonian limit .
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
double * p_ray_eq2_pi
Coordinate radius at , .
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
virtual const Itbl & l_surf() const
Description of the surface of fluid 1: returns a 2-D Itbl containing the values of the domain index ...
virtual double mom_quad_old() const
Part of the quadrupole moment.
virtual double grv2() const
Error on the virial identity GRV2.
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
double * p_mass_g
Gravitational mass.
Definition: etoile.h:551
virtual double mean_radius2() const
Mean radius for fluid 2.
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
double ray_pole2() const
Coordinate radius for fluid 2 at [r_unit].
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
virtual double angu_mom_1() const
Angular momentum of fluid 1.
Tenseur a_car
Total conformal factor .
Definition: etoile.h:518
double * p_angu_mom_2
Angular momentum of fluid 2.
Tenseur ent2
Log-enthalpy for the second fluid.
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
Tenseur sphph_euler
The component of the stress tensor .
virtual double mom_quad() const
Quadrupole moment.
int get_dzpuis() const
Returns dzpuis.
Definition: cmp.h:903
double * p_ray_eq2_pis2
Coordinate radius at , .
double * p_ray_pole2
Coordinate radius at .
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:460
virtual double r_circ2() const
Circumferential radius for fluid 2.
Tenseur nbar2
Baryon density in the fluid frame, for fluid 2.
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
virtual double coupling_mominert_1() const
Quantities used to study the different fluid couplings: , and .
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 angu_mom() const
Angular momentum.
double mass_b2() const
Baryon mass of fluid 2.
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 alpha_eos
Coefficient relative to entrainment effects.
Tenseur delta_car
The "relative velocity" (squared) of the two fluids.
double ray_eq2() const
Coordinate radius for fluid 2 at , [r_unit].
Tenseur & dzeta
Metric potential = beta_auto.
Definition: etoile.h:1537
virtual double aplat2() const
Flatening r_pole/r_eq for fluid 2.
double * p_aplat2
Flatening r_pole/r_eq of fluid no.2.
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:464
double * p_mass_b
Baryon mass.
Definition: etoile.h:550
double get_m1() const
Return the individual particule mass .
Definition: eos_bifluid.h:263
virtual double area2() const
Surface area for fluid 2.
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
double * p_angu_mom
Angular momentum.
Definition: etoile.h:1634
Itbl * p_l_surf2
Description of the surface of fluid 2: 2-D Itbl containing the values of the domain index l on the su...
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