LORENE
hot_strot_cfc_equil.C
1 /*
2  * Function Hot_star_rot_CFC::equilibrium
3  *
4  * (see file hot_star_rot_cfc.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2024 Jerome Novak
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License version 2
15  * as published by the Free Software Foundation.
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 // Lorene headers
30 #include "hot_star_rot_cfc.h"
31 
32 #include "utilitaires.h"
33 #include "unites.h"
34 
35 namespace Lorene {
36 void Hot_star_rot_CFC::equilibrium(double ent_c, double omega0,
37  double fact_omega, int , const Tbl& ent_limit,
38  const Itbl& icontrol, const Tbl& control,
39  double mbar_wanted, double aexp_mass, Tbl& diff){
40 
41 
42  // Fundamental constants and units
43  // --------------------------------
44  using namespace Unites ;
45 
46  // For the display
47  // ---------------
48  char display_bold[]="x[1m" ; display_bold[0] = 27 ;
49  char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
50 
51  // Grid parameters
52  // ----------------
53 
54  const Mg3d* mg = mp.get_mg() ;
55  int nz = mg->get_nzone() ; // total number of domains
56  int nzm1 = nz - 1 ;
57 
58  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
59  int type_t = mg->get_type_t() ;
60  assert( ( type_t == SYM) || (type_t == NONSYM) ) ;
61  int l_b = nzet - 1 ;
62  int i_b = mg->get_nr(l_b) - 1 ;
63  int j_b = (type_t == SYM ? mg->get_nt(l_b) - 1 : mg->get_nt(l_b)/2 ) ;
64  int k_b = 0 ;
65 
66  // Value of the enthalpy defining the surface of the star
67  double ent_b = ent_limit(nzet-1) ;
68 
69  // Parameters to control the iteration
70  // -----------------------------------
71 
72  int mer_max = icontrol(0) ;
73  int mer_rot = icontrol(1) ;
74  int mer_change_omega = icontrol(2) ;
75  int mer_fix_omega = icontrol(3) ;
76  int mer_mass = icontrol(4) ;
77  int delta_mer_kep = icontrol(5) ;
78 
79  // Protections:
80  if (mer_change_omega < mer_rot) {
81  cout << "Hot_star_rot_CFC::equilibrium: mer_change_omega < mer_rot !"
82  << endl ;
83  cout << " mer_change_omega = " << mer_change_omega << endl ;
84  cout << " mer_rot = " << mer_rot << endl ;
85  abort() ;
86  }
87  if (mer_fix_omega < mer_change_omega) {
88  cout << "Hot_star_rot_CFC::equilibrium: mer_fix_omega < mer_change_omega !"
89  << endl ;
90  cout << " mer_fix_omega = " << mer_fix_omega << endl ;
91  cout << " mer_change_omega = " << mer_change_omega << endl ;
92  abort() ;
93  }
94 
95  // In order to converge to a given baryon mass, shall the central
96  // enthalpy be varied or Omega ?
97  bool change_ent = true ;
98  if (mer_mass < 0) {
99  change_ent = false ;
100  mer_mass = abs(mer_mass) ;
101  }
102 
103 
104  double precis = control(0) ;
105  double omega_ini = control(1) ;
106  double relax = control(2) ;
107  double relax_prev = double(1) - relax ;
108 
109  // Error indicators
110  // ----------------
111 
112  diff.annule_hard() ;
113  double& diff_ent = diff.set(0) ;
114 
115  double alpha_r = 1 ;
116 
117  // Initializations
118  // ---------------
119 
120  // Initial angular velocities
121  omega = 0 ;
122 
123  double accrois_omega = (omega0 - omega_ini) /
124  double(mer_fix_omega - mer_change_omega) ;
125 
126 
127  update_metric() ; //update of the metric quantities
128 
129  update_entropy() ; // update entropy
130 
131  equation_of_state() ; // update of the density, pressure,...etc
132 
133  hydro_euler() ; //update of the hydro quantities relative to the
134  // Eulerian observer
135 
136  // Quantities at the previous step :
137  Scalar ent_prev = ent ;
138  Scalar logn_prev = logn ;
139  Scalar psi_prev = psi ;
140  // Vector beta_prev = beta ;
141 
142  // Output files
143  // -------------
144 
145  ofstream fichconv("convergence.d") ; // Output file for diff_ent
146  fichconv << "# diff_ent GRV2 max_triax vit_triax" << endl ;
147 
148  ofstream fichfreq("frequency.d") ; // Output file for omega
149  fichfreq << "# f [Hz]" << endl ;
150 
151  ofstream fichevol("evolution.d") ; // Output file for various quantities
152  fichevol <<
153  "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
154  << endl ;
155 
156  diff_ent = 1 ;
157  double err_grv2 = 1 ;
158 
159 
160 
161 //=========================================================================
162 // Start of iteration
163 //=========================================================================
164 
165  for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++) {
166 
167  cout << "-----------------------------------------------" << endl ;
168  cout << "step: " << mer << endl ;
169  cout << "diff_ent = " << display_bold << diff_ent << display_normal
170  << endl ;
171  cout << "err_grv2 = " << err_grv2 << endl ;
172  fichconv << mer ;
173  fichfreq << mer ;
174  fichevol << mer ;
175 
176 
177  // switch on rotation
178  if (mer >= mer_rot) {
179 
180  if (mer < mer_change_omega) {
181  omega = omega_ini ;
182  }
183  else {
184  if (mer <= mer_fix_omega) {
185  omega = omega_ini + accrois_omega *
186  (mer - mer_change_omega) ;
187  }
188  }
189 
190 
191  }
192 
193  //---------------------------------------------//
194  // Resolution of the Poisson equation for psi //
195  //---------------------------------------------//
196 
197  Scalar psi_new(mp) ;
198 
199  solve_psi( psi_new ) ;
200 
201  psi_new.std_spectral_base() ;
202 
203 
204  //---------------------------------------------------//
205  // Resolution of the Poisson equation for logn //
206  // Note: ln_f is due to the fluid part //
207  // ln_q is due to the quadratic metric part //
208  //---------------------------------------------------//
209 
210  Scalar ln_f_new(mp) ;
211  Scalar ln_q_new(mp) ;
212 
213  solve_logn_f( ln_f_new ) ;
214  solve_logn_q( ln_q_new ) ;
215 
216  ln_f_new.std_spectral_base() ;
217  ln_q_new.std_spectral_base() ;
218 
219 
220  //--------------------------------------------------//
221  // Resolution of the Poisson equation for shift //
222  //--------------------------------------------------//
223 
224  Vector beta_new(mp, CON, mp.get_bvect_spher()) ;
225 
226  solve_shift( beta_new ) ;
227 
228  //------------------------------------
229  // Determination of the fluid velocity
230  //------------------------------------
231 
232  if (mer > mer_fix_omega + delta_mer_kep) {
233 
234  omega *= fact_omega ; // Increase of the angular velocity if
235  } // fact_omega != 1
236 
237  bool omega_trop_grand = false ;
238  bool kepler = true ;
239 
240  while ( kepler ) {
241 
242  // Possible decrease of Omega to ensure a velocity < c
243 
244  bool superlum = true ;
245 
246  while ( superlum ){
247 
248  // New fluid velocity :
249  //
250 
251  u_euler.set(1).set_etat_zero() ;
252  u_euler.set(2).set_etat_zero() ;
253 
254  u_euler.set(3) = omega ;
255  u_euler.set(3).annule(nzet,nzm1) ; // nzet is defined in class Hot_star
257  u_euler.set(3).mult_rsint() ;
258  u_euler.set(3) += beta(3) ;
259  u_euler.set(3).annule(nzet,nzm1) ;
260 
261  u_euler = u_euler / nn ;
262 
263 
264  // v2 (square of norm of u_euler)
265  // -------------------------------
266 
267  v2 = contract(contract(gamma.cov(), 0, u_euler, 0), 0, u_euler, 0) ;
268 
269  // Is the new velocity larger than c in the equatorial plane ?
270 
271  superlum = false ;
272 
273  for (int l=0; l<nzet; l++) {
274  for (int i=0; i<mg->get_nr(l); i++) {
275 
276  double u2 = v2.val_grid_point(l, 0, j_b, i) ;
277  if (u2 >= 1.) { // superluminal velocity
278  superlum = true ;
279  cout << "U > c for l, i : " << l << " " << i
280  << " U = " << sqrt( u2 ) << endl ;
281  }
282  }
283  }
284  if ( superlum ) {
285  cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ;
286  omega /= fact_omega ; // Decrease of Omega
287  cout << "New rotation frequency : "
288  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
289  omega_trop_grand = true ;
290  }
291  } // end of while ( superlum )
292 
293 
294  // New computation of U (this time is not superluminal)
295  // as well as of gam_euler, ener_euler,...etc
296 
297  hydro_euler() ;
298 
299 
300 
301  //--------------------------------//
302  // First integral of motion //
303  //--------------------------------//
304 
305  Scalar mlngamma(mp) ; // -log( gam_euler )
306 
307  mlngamma = - log( gam_euler ) ;
308 
309  update_entropy() ;
310  Scalar expmHT(mp) ;
311  expmHT = exp(-ent)*temp ;
312  expmHT.std_spectral_base() ;
313  Scalar prim_entropy = expmHT ;
314  prim_entropy *= entropy.dsdr() ;
315  prim_entropy.annule(nzet, nzm1) ;
316 
317  prim_entropy = prim_entropy.primr(false) ;
318 
319  // Equatorial values of various potentials :
320  double ln_f_b = ln_f_new.val_grid_point(l_b, k_b, j_b, i_b) ;
321  double ln_q_b = ln_q_new.val_grid_point(l_b, k_b, j_b, i_b) ;
322  double mlngamma_b = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ;
323  double prim_entr_b = prim_entropy.val_grid_point(l_b, k_b, j_b, i_b) ;
324 
325  // Central values of various potentials :
326  double ln_f_c = ln_f_new.val_grid_point(0,0,0,0) ;
327  double ln_q_c = ln_q_new.val_grid_point(0,0,0,0) ;
328  double mlngamma_c = 0 ;
329  double prim_entr_c = prim_entropy.val_grid_point(0,0,0,0) ;
330 
331  // Scale factor to ensure that the (log of) enthalpy is equal to
332  // ent_b at the equator
333  double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
334  - prim_entr_c + prim_entr_b + ln_q_c - ln_q_b) / ( ln_f_b - ln_f_c ) ;
335 
336  alpha_r = sqrt(alpha_r2) ;
337 
338  cout << "alpha_r = " << alpha_r << endl ;
339 
340  // Rescaling of the grid (no adaptation!)
341  //---------------------------------------
342  mp.homothetie(alpha_r) ;
343 
344  // Readjustment of logn :
345  // -----------------------
346 
347  logn = alpha_r2 * ln_f_new + ln_q_new ;
348 
349  double logn_c = logn.val_grid_point(0,0,0,0) ;
350 
351  // First integral of motion -> (log of) enthalpy in all space
352  // ----------------------------------------------------------
353 
354  ent = (ent_c + logn_c + mlngamma_c - prim_entr_c) - logn - mlngamma + prim_entropy ;
355 
356  //---------------------------------
357  // Equation of state
358  //---------------------------------
359 
360  equation_of_state() ; // computes new values for nbar (n), ener (e),
361  // and press (p) from the new ent (H)
362 
363  // --------------------------------------------------------------
364  // Test: is the enthalpy negative somewhere in the equatorial plane
365  // inside the star?
366  // --------------------------------------------------------
367 
368  kepler = false ;
369  for (int l=0; l<nzet; l++) {
370  int imax = mg->get_nr(l) - 1 ;
371  if (l == l_b) imax-- ; // The surface point is skipped
372  for (int i=0; i<imax; i++) {
373  if ( ent.val_grid_point(l, 0, j_b, i) < 0. ) {
374  kepler = true ;
375  cout << "ent < 0 for l, i : " << l << " " << i
376  << " ent = " << ent.val_grid_point(l, 0, j_b, i) << endl ;
377  }
378  }
379  }
380 
381  if ( kepler ) {
382  cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
383  omega /= fact_omega ; // Omega is decreased
384  cout << "New rotation frequency : "
385  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
386  omega_trop_grand = true ;
387  }
388 
389  } // End of while ( kepler )
390 
391  if ( omega_trop_grand ) { // fact_omega is decreased for the
392  // next step
393  fact_omega = sqrt( fact_omega ) ;
394  cout << "**** New fact_omega : " << fact_omega << endl ;
395  }
396 
397  hydro_euler() ;
398 
399  beta = beta_new ;
400 
401  //---------------------------------------
402  // Calculate error of the GRV2 identity
403  //---------------------------------------
404 
405  err_grv2 = grv2() ;
406 
407 
408  //--------------------------------------
409  // Relaxations on some quantities....?
410  //
411  //--------------------------------------
412 
413  if (mer >= 10) {
414  logn = relax * logn + relax_prev * logn_prev ;
415 
416  psi = relax * psi_new + relax_prev * psi_prev ;
417 
418  }
419 
420  // Update of the metric quantities :
421 
422  update_metric() ;
423 
424  //---------------------------
425  // Informations display
426  // More to come later......
427  //---------------------------
428 
429  // partial_display(cout) ; // What is partial_display(cout) ?
430  fichfreq << " " << omega / (2*M_PI) * f_unit ;
431  fichevol << " " << ent_c ;
432 
433 
434  //-----------------------------------------
435  // Convergence towards a given baryon mass
436  //-----------------------------------------
437 
438  if (mer > mer_mass) {
439 
440  double xx ;
441  if (mbar_wanted > 0.) {
442  xx = mass_b() / mbar_wanted - 1. ;
443  cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx
444  << endl ;
445  }
446  else{
447  xx = mass_g() / fabs(mbar_wanted) - 1. ;
448  cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx
449  << endl ;
450  }
451  double xprog = ( mer > 2*mer_mass) ? 1. :
452  double(mer-mer_mass)/double(mer_mass) ;
453  xx *= xprog ;
454  double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
455  double fact = pow(ax, aexp_mass) ;
456  cout << " xprog, xx, ax, fact : " << xprog << " " <<
457  xx << " " << ax << " " << fact << endl ;
458 
459  if ( change_ent ) {
460  ent_c *= fact ;
461  }
462  else {
463  if (mer%4 == 0) omega *= fact ;
464  }
465  }
466 
467 
468  //-----------------------------------------------------------
469  // Relative change in enthalpy with respect to previous step
470  // ** Check: Is diffrel(ent, ent_prev) ok?
471  //-----------------------------------------------------------
472 
473  Tbl diff_ent_tbl = diffrel( ent, ent_prev ) ;
474  diff_ent = diff_ent_tbl(0) ;
475  for (int l=1; l<nzet; l++) {
476  diff_ent += diff_ent_tbl(l) ;
477  }
478  diff_ent /= nzet ;
479 
480  fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
481  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
482 
483  //------------------------------
484  // Recycling for the next step
485  //------------------------------
486 
487  ent_prev = ent ;
488  logn_prev = logn ;
489  psi_prev = psi ;
490 
491  fichconv << endl ;
492  fichfreq << endl ;
493  fichevol << endl ;
494  fichconv.flush() ;
495  fichfreq.flush() ;
496  fichevol.flush() ;
497 
498 
499  } // End of main loop
500 
501  //=================================================
502  // End of iteration
503  //=================================================
504 
505  fichconv.close() ;
506  fichfreq.close() ;
507  fichevol.close() ;
508 
509 
510 }
511 }
void solve_shift(Vector &shift_new) const
Solution of the shift equation for rotating stars in CFC.
Metric gamma
3-metric
Definition: hot_star.h:140
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
Scalar temp
Temperature field (in MeV)
Definition: hot_star.h:98
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
Scalar psi
Conformal factor .
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
Definition: scalar.C:397
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:809
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
virtual double grv2() const
Error on the virial identity GRV2.
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:402
void solve_logn_q(Scalar &ln_q_new) const
Solution of the quadratic part of the Poisson equation for the lapse for rotating stars in CFC...
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: hot_star.C:489
Scalar primr(bool null_infty=true) const
Computes the radial primitive which vanishes for .
Definition: scalar_integ.C:106
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition: grilles.h:504
Basic integer array class.
Definition: itbl.h:122
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
Scalar entropy
Entropy per baryon field (in $k_B$)
Definition: hot_star.h:99
Scalar nn
Lapse function N .
Definition: hot_star.h:130
Tensor field of valence 1.
Definition: vector.h:188
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
void solve_logn_f(Scalar &ln_f_new) const
Solution of the `‘matter’&#39; part of the Poisson equation for the lapse for rotating stars in CFC...
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:652
Scalar logn
Logarithm of the lapse N .
Definition: hot_star.h:127
int nzet
Number of domains of *mp occupied by the star.
Definition: hot_star.h:84
virtual double mass_b() const
Baryonic mass.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
virtual void homothetie(double lambda)=0
Sets a new radial scale.
double omega
Rotation angular velocity ([f_unit] )
void equilibrium(double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff)
Computes an equilibrium configuration.
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: hot_star.h:109
virtual double mass_g() const
Gravitational mass.
Scalar ent
Log-enthalpy.
Definition: hot_star.h:91
Map & mp
Mapping associated with the star.
Definition: hot_star.h:81
void solve_psi(Scalar &psi_new)
Solution of the equations for the conformal factor for rotating stars in CFC.
void mult_rsint()
Multiplication by everywhere; dzpuis is not changed.
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 .
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
Definition: metric.C:283
Multi-domain grid.
Definition: grilles.h:276
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition: cmp_math.C:325
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
void update_metric()
Computes metric quantities from known potentials.
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:112
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:476
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
void update_entropy()
Updates the part of the entropy from entropy_param.
void annule_hard()
Sets the Tbl to zero in a hard way.
Definition: tbl.C:375
Vector beta
Shift vector.
Definition: hot_star.h:133
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: hot_star.h:112