LORENE
hot_strot_diff_cfc_equil.C
1 /*
2  * Function Hot_star_rot_diff_CFC::equilibrium
3  *
4  * (see file hot_star_rot_diff_cfc.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2025 Santiago Jaraba
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 // C headers
29 #include <cmath>
30 
31 // Lorene headers
32 #include "hot_star_rot_diff_cfc.h"
33 #include "param.h"
34 #include "graphique.h"
35 #include "utilitaires.h"
36 #include "unites.h"
37 
38 
39 namespace Lorene{
40 void Hot_star_rot_diff_CFC::equilibrium(double ent_c, double omega_c0, double fact_omega,
41  int , const Tbl& ent_limit,
42  const Itbl& icontrol,
43  const Tbl& control, double mbar_wanted,
44  double aexp_mass, Tbl& diff, Param*){
45 
46  // Fundamental constants and units
47  // --------------------------------
48  using namespace Unites ;
49 
50  // For the display
51  // ---------------
52  char display_bold[]="x[1m" ; display_bold[0] = 27 ;
53  char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
54 
55 
56  // Grid parameters
57  // ----------------
58 
59  const Mg3d* mg = mp.get_mg() ;
60  int nz = mg->get_nzone() ; // total number of domains
61  int nzm1 = nz - 1 ;
62 
63  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
64  int type_t = mg->get_type_t() ;
65  assert( ( type_t == SYM) || (type_t == NONSYM) ) ;
66  int l_b = nzet - 1 ;
67  int i_b = mg->get_nr(l_b) - 1 ;
68  int j_b = (type_t == SYM ? mg->get_nt(l_b) - 1 : mg->get_nt(l_b)/2 ) ;
69  int k_b = 0 ;
70 
71  // Value of the enthalpy defining the surface of the star
72  double ent_b = ent_limit(nzet-1) ;
73 
74  // Parameters to control the iteration
75  // -----------------------------------
76 
77  int mer_max = icontrol(0) ;
78  int mer_rot = icontrol(1) ;
79  int mer_change_omega = icontrol(2) ;
80  int mer_fix_omega = icontrol(3) ;
81  int mer_mass = icontrol(4) ;
82  int mer_triax = icontrol(5) ;
83  int delta_mer_kep = icontrol(6) ;
84 
85  // Protections:
86  if (mer_change_omega < mer_rot) {
87  cout << "Hot_star_rot_diff_CFC::equilibrium: mer_change_omega < mer_rot !" << endl ;
88  cout << " mer_change_omega = " << mer_change_omega << endl ;
89  cout << " mer_rot = " << mer_rot << endl ;
90  abort() ;
91  }
92  if (mer_fix_omega < mer_change_omega) {
93  cout << "Hot_star_rot_diff_CFC::equilibrium: mer_fix_omega < mer_change_omega !"
94  << endl ;
95  cout << " mer_fix_omega = " << mer_fix_omega << endl ;
96  cout << " mer_change_omega = " << mer_change_omega << endl ;
97  abort() ;
98  }
99 
100  // In order to converge to a given baryon mass, shall the central
101  // enthalpy be varied or Omega ?
102  bool change_ent = true ;
103  if (mer_mass < 0) {
104  change_ent = false ;
105  mer_mass = abs(mer_mass) ;
106  }
107 
108  double precis = control(0) ;
109  double omega_ini = control(1) ;
110  double relax = control(2) ;
111  double relax_prev = double(1) - relax ;
112  double ampli_triax = control(3) ;
113 
114 
115  // Error indicators
116  // ----------------
117 
118  diff.annule_hard() ;
119  double& diff_ent = diff.set(0) ;
120  double& vit_triax = diff.set(7) ;
121 
122  double alpha_r = 1 ;
123 
124  // Initializations
125  // ---------------
126 
127  // Initial central angular velocity
128  double omega_c = 0 ;
129 
130  double accrois_omega = (omega_c0 - omega_ini) /
131  double(mer_fix_omega - mer_change_omega) ;
132 
133 
134  update_metric() ; //update of the metric quantities
135 
136  update_entropy() ; // update of the entropy
137 
138  equation_of_state() ; // update of the density, pressure,...etc
139 
140  hydro_euler() ; //update of the hydro quantities relative to the
141  // Eulerian observer
142 
143  // Quantities at the previous step :
144  Scalar ent_prev = ent ;
145  Scalar logn_prev = logn ;
146  Scalar psi_prev = psi ;
147  Vector beta_prev = beta ;
148 
149  // Output files
150  // -------------
151 
152  ofstream fichconv("convergence.d") ; // Output file for diff_ent
153  fichconv << "# diff_ent GRV2 max_triax vit_triax" << endl ;
154 
155  ofstream fichfreq("frequency.d") ; // Output file for omega_c
156  fichfreq << "# f [Hz]" << endl ;
157 
158  ofstream fichevol("evolution.d") ; // Output file for various quantities
159  fichevol <<
160  "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
161  << endl ;
162 
163  diff_ent = 1 ;
164  double err_grv2 = 1 ;
165  double max_triax_prev = 0 ; // Triaxial amplitude at previous step
166 
167 
168  //=========================================================================
169  // Start of iteration
170  //=========================================================================
171 
172  for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
173 
174  cout << "-----------------------------------------------" << endl ;
175  cout << "step: " << mer << endl ;
176  cout << "diff_ent = " << display_bold << diff_ent << display_normal
177  << endl ;
178  cout << "err_grv2 = " << err_grv2 << endl ;
179  fichconv << mer ;
180  fichfreq << mer ;
181  fichevol << mer ;
182 
183  // Switch on rotation
184  if (mer >= mer_rot) {
185 
186  if (mer < mer_change_omega) {
187  omega_c = omega_ini ;
188  }
189  else {
190  if (mer <= mer_fix_omega) {
191  omega_c = omega_ini + accrois_omega *
192  (mer - mer_change_omega) ;
193  }
194  }
195 
196 
197  }
198 
199  //---------------------------------------------//
200  // Resolution of the Poisson equation for psi //
201  //---------------------------------------------//
202 
203  Scalar psi_new(mp) ;
204 
205  solve_psi( psi_new ) ;
206 
207  psi_new.std_spectral_base() ;
208 
209 
210  //---------------------------------------------------//
211  // Resolution of the Poisson equation for logn //
212  // Note: ln_f is due to the fluid part //
213  // ln_q is due to the quadratic metric part //
214  //---------------------------------------------------//
215 
216  Scalar ln_f_new(mp) ;
217  Scalar ln_q_new(mp) ;
218 
219  solve_logn_f( ln_f_new ) ;
220  solve_logn_q( ln_q_new ) ;
221 
222  ln_f_new.std_spectral_base() ;
223  ln_q_new.std_spectral_base() ;
224 
225  //---------------------------------------
226  // Triaxial perturbation of ln_ff
227  //---------------------------------------
228 
229  if (mer == mer_triax) {
230 
231  if ( mg->get_np(0) == 1 ) {
232  cout <<
233  "Hot_star_rot_diff::equilibrium: np must be stricly greater than 1"
234  << endl << " to set a triaxial perturbation !" << endl ;
235  abort() ;
236  }
237 
238  const Coord& phi = mp.phi ;
239  const Coord& sint = mp.sint ;
240  Scalar perturb(mp) ;
241  perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
242  ln_f_new = ln_f_new * perturb ;
243 
244  ln_f_new.std_spectral_base() ; // set the bases for spectral expansions
245  // to be the standard ones for a
246  // scalar field
247 
248  }
249 
250  // Monitoring of the triaxial perturbation
251  // ---------------------------------------
252 
253  const Valeur& va_nuf = ln_f_new.get_spectral_va() ;
254  va_nuf.coef() ; // Computes the spectral coefficients
255  double max_triax = 0 ;
256 
257  if ( mg->get_np(0) > 1 ) {
258 
259  for (int l=0; l<nz; l++) { // loop on the domains
260  for (int j=0; j<mg->get_nt(l); j++) {
261  for (int i=0; i<mg->get_nr(l); i++) {
262 
263  // Coefficient of cos(2 phi) :
264  double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ;
265 
266  // Coefficient of sin(2 phi) :
267  double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ;
268 
269  double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
270 
271  max_triax = ( xx > max_triax ) ? xx : max_triax ;
272  }
273  }
274  }
275 
276  }
277 
278  cout << "Triaxial part of nuf : " << max_triax << endl ;
279 
280 
281  //--------------------------------------------------//
282  // Resolution of the Poisson equation for shift //
283  //--------------------------------------------------//
284 
285  Vector beta_new(mp, CON, mp.get_bvect_spher()) ;
286 
287  solve_shift( beta_new ) ;
288 
289  //-----------------------------------------
290  // Determination of the fluid velociy U
291  //-----------------------------------------
292 
293  if (mer > mer_fix_omega + delta_mer_kep) {
294 
295  omega_c *= fact_omega ; // Increase of the angular velocity if
296  } // fact_omega != 1
297 
298 
299  bool omega_trop_grand = false ;
300  bool kepler = true ;
301 
302  while ( kepler ) {
303 
304  // Possible decrease of Omega to ensure a velocity < c
305 
306  bool superlum = true ;
307 
308  while ( superlum ) {
309 
310  // Computation of Omega(r,theta)
311 
312  if (omega_c == 0.) {
313  omega_field = 0 ;
314  }
315  else if (par_frot.get_taille() == 5) { // Uryu rotation profile
316  if (par_frot(1) == 1. && par_frot(2) == 1.) { // Rigid rotation limit in Uryu profile
317  omega_field = omega_c ;
318  for (int l=nzet+1; l<nz; l++) { // omega_field set to 0 far from the star
319  omega_field.set_domain(l) = 0 ;
320  }
321  prim_field = 0 ;
322  omega_min = omega_c ;
323  omega_max = omega_c ;
324  }
325  else { // Usual case
326  double lambda1 = par_frot(1) ;
327  double lambda2 = par_frot(2) ;
328  int p = 1 ;
329  int q = 3 ;
330 
331  // First, compute F_e and F_max
332  Scalar rsint(mp) ;
333  rsint = 1 ;
334  rsint.annule_domain(nzm1) ;
335  rsint.std_spectral_base() ;
336  rsint.mult_rsint() ;
337  Scalar nnn2 = nn * nn ;
338 
339  double omnp = omega_field.val_grid_point(l_b, k_b, j_b, i_b)*rsint.val_grid_point(l_b, k_b, j_b, i_b)
340  + beta(3).val_grid_point(l_b, k_b, j_b, i_b) ;
341  double F_e = psi4.val_grid_point(l_b, k_b, j_b, i_b)*rsint.val_grid_point(l_b, k_b, j_b, i_b) * omnp /
342  ( nnn2.val_grid_point(l_b, k_b, j_b, i_b) -
343  psi4.val_grid_point(l_b, k_b, j_b, i_b) * omnp * omnp ) ;
344 
345  if (F_e != double(0)) { // Not the first iteration
346  int l, k, j, i ;
347  double om_max = 0 ;
348  double om_lkji = 0 ;
349  double F_max = 0 ;
350  for (l=0; l<nzet+1; l++) {
351  for (k=0; k<mg->get_np(l); k++) {
352  for (j=0; j<mg->get_nt(l); j++) {
353  for (i=0; i<mg->get_nr(l); i++) {
354  om_lkji = omega_field.val_grid_point(l, k, j, i) ;
355  if (om_lkji > om_max) {
356  om_max = om_lkji ;
357  omnp = om_max*rsint.val_grid_point(l, k, j, i) + beta(3).val_grid_point(l, k, j, i) ;
358  F_max = psi4.val_grid_point(l, k, j, i)*rsint.val_grid_point(l, k, j, i) * omnp /
359  ( nnn2.val_grid_point(l, k, j, i) -
360  psi4.val_grid_point(l, k, j, i) * omnp * omnp ) ;
361  }
362  }
363  }
364  }
365  }
366 
367  // Now A^2 and B^2, parameters for omega_frot, can be computed
368 
369  if (lambda2 == -1.) { // Special flag to keep A^2 and B^2 constant
370  // Computing lambda1 from A^2 and B^2
371  double A2 = par_frot(3) ;
372  double B2 = par_frot(4) ;
373  lambda1 = (1 + pow(F_max/(B2*omega_c), p))/(1 + pow(F_max/(A2*omega_c), q+p)) ;
374  }
375  else if (F_e < 0) { // Work in progress: fixing issues with A^2 and B^2 parameters
376  cout << "Warning: F_e=" << F_e << " < 0. Aborting." << endl ;
377  abort() ;
378  //F_e = 1.001*pow(lambda1/lambda2, 1./q)*F_max ;
379  //cout << F_e << endl ;
380  par_frot.set(3) = 10000. ;
381  par_frot.set(4) = 10000. ;
382  }
383  /*else if (F_e < F_max) {
384  cout << "Warning: F_e=" << F_e << " < F_max=" << F_max << ". Impossible to fix by changing lambda1 and lambda2." << endl ;
385  cout << "Trying with the parameters of the previous iteration" << endl ;
386  }*/
387  else {
388  /*while (F_e < pow(lambda1/lambda2, 1./q)*F_max) {
389  // This would lead to unphysical A^2 and B^2. Relaxing lambda1 and lambda2 only in this iteration
390  lambda1 *= .9 ;
391  lambda2 *= 1.1 ;
392  // lambda2 = 1.1*lambda1*pow(F_max/F_e, q) ;
393  cout << "Warning: unphysical A^2 and B^2. Relaxing lambda1=" << par_frot(1)
394  << "->" << lambda1 << ", lambda2=" << par_frot(2) << "->" << lambda2 << endl ;
395 
396  }*/
397  //cout << lambda2*pow(F_e, q)-lambda1*pow(F_max, q) << endl ;
398  par_frot.set(3) = pow(pow(F_e*F_max,p)/pow(omega_c, p+q) * (lambda2*pow(F_e, q)-lambda1*pow(F_max, q) )
399  / ((lambda1-1)*pow(F_e, p) - (lambda2-1)*pow(F_max, p)), 1./(p+q) ) ;
400  par_frot.set(4) = F_e*F_max/omega_c * pow((lambda2*pow(F_e, q)-lambda1*pow(F_max, q))
401  / (lambda2*(lambda1-1)*pow(F_e, p+q) - lambda1*(lambda2-1)*pow(F_max, p+q)), 1./p) ;
402 
403  cout << F_e << " " << F_max << endl ;
404  }
405 
406  }
407 
408  par_frot.set(0) = omega_c ;
409 
410  cout << "Parameters of Omega(F) : " << par_frot << endl ;
411 
412  double omeg_min = lambda2 * omega_c ;
413  double omeg_max = lambda1 * omega_c ;
414  double precis1 = 1.e-14 ;
415  int nitermax1 = 1000 ;
416 
417  fait_omega_field(omeg_min, omeg_max, precis1, nitermax1) ;
418  }
419  }
420  else { // Usual rotation profile
421  par_frot.set(0) = omega_c ;
422  if (par_frot(2) != double(0)) { // fixed a = R_eq / R_0
423  par_frot.set(1) = ray_eq() / par_frot(2) ;
424  }
425  double omeg_min = 0 ;
426  double omeg_max = omega_c ;
427  double precis1 = 1.e-14 ;
428  int nitermax1 = 100 ;
429 
430  fait_omega_field(omeg_min, omeg_max, precis1, nitermax1) ;
431  }
432 
433  // New fluid velocity :
434  //
435 
436  u_euler.set(1).set_etat_zero() ;
437  u_euler.set(2).set_etat_zero() ;
438 
439  u_euler.set(3) = omega_field ;
440  u_euler.set(3).annule(nzet,nzm1) ; // nzet is defined in class Hot_star
442  u_euler.set(3).mult_rsint() ;
443  u_euler.set(3) += beta(3) ;
444  u_euler.set(3).annule(nzet,nzm1) ;
445 
446  u_euler = u_euler / nn ;
447 
448 
449  // v2 (square of norm of u_euler)
450  // -------------------------------
451 
452  v2 = contract(contract(gamma.cov(), 0, u_euler, 0), 0, u_euler, 0) ;
453 
454  // Is the new velocity larger than c in the equatorial plane ?
455 
456  superlum = false ;
457 
458  for (int l=0; l<nzet; l++) {
459  for (int i=0; i<mg->get_nr(l); i++) {
460 
461  double u2 = v2.val_grid_point(l, 0, j_b, i) ;
462  if (u2 >= 1.) { // superluminal velocity
463  superlum = true ;
464  cout << "U > c for l, i : " << l << " " << i
465  << " U = " << sqrt(u2) << endl ;
466  }
467  }
468  }
469  if ( superlum ) {
470  cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ;
471  omega_c /= fact_omega ; // Decrease of Omega_c
472  cout << "New central rotation frequency : "
473  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
474  omega_trop_grand = true ;
475  }
476  } // end of while ( superlum )
477 
478 
479  // New computation of U (which this time is not superluminal)
480  // as well as of gam_euler, ener_euler, etc...
481  // -----------------------------------
482 
483  hydro_euler() ;
484 
485 
486  //------------------------------------------------------
487  // First integral of motion
488  //------------------------------------------------------
489 
490  Scalar mlngamma(mp) ; // -log( gam_euler )
491 
492  mlngamma = - log( gam_euler ) ;
493 
494  update_entropy() ;
495  Scalar expmHT(mp) ;
496  expmHT = exp(-ent)*temp ;
497  expmHT.std_spectral_base() ;
498  Scalar prim_entropy = expmHT ;
499  prim_entropy *= entropy.dsdr() ;
500  prim_entropy.annule(nzet, nzm1) ;
501 
502  prim_entropy = prim_entropy.primr(false) ;
503 
504  // Equatorial values of various potentials :
505  double ln_f_b = ln_f_new.val_grid_point(l_b, k_b, j_b, i_b) ;
506  double ln_q_b = ln_q_new.val_grid_point(l_b, k_b, j_b, i_b) ;
507  double mlngamma_b = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ;
508  double prim_entr_b = prim_entropy.val_grid_point(l_b, k_b, j_b, i_b) ;
509  double primf_b = prim_field.val_grid_point(l_b, k_b, j_b, i_b) ;
510 
511 
512  // Central values of various potentials :
513  double ln_f_c = ln_f_new.val_grid_point(0,0,0,0) ;
514  double ln_q_c = ln_q_new.val_grid_point(0,0,0,0) ;
515  double mlngamma_c = 0 ;
516  double prim_entr_c = prim_entropy.val_grid_point(0,0,0,0) ;
517  double primf_c = prim_field.val_grid_point(0,0,0,0) ;
518 
519  // Scale factor to ensure that the (log of) enthalpy is equal to
520  // ent_b at the equator
521  double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
522  - prim_entr_c + prim_entr_b + ln_q_c - ln_q_b + primf_c - primf_b)
523  / ( ln_f_b - ln_f_c ) ;
524  alpha_r = sqrt(alpha_r2) ;
525 
526  cout << "alpha_r = " << alpha_r << endl ;
527 
528  // Rescaling of the grid (no adaptation!)
529  //---------------------------------------
530  mp.homothetie(alpha_r) ;
531 
532  // Readjustment of logn :
533  // -----------------------
534 
535  logn = alpha_r2 * ln_f_new + ln_q_new ;
536 
537  double logn_c = logn.val_grid_point(0,0,0,0) ;
538 
539  // First integral of motion -> (log of) enthalpy in all space
540  // ----------------------------------------------------------
541 
542  ent = (ent_c + logn_c + mlngamma_c - prim_entr_c)
543  - logn - mlngamma - prim_field + prim_entropy ;
544 
545  //---------------------------------
546  // Equation of state
547  //---------------------------------
548 
549  equation_of_state() ; // computes new values for nbar (n), ener (e),
550  // and press (p) from the new ent (H)
551 
552  // Test: is the enthalpy negative somewhere in the equatorial plane
553  // inside the star?
554  // --------------------------------------------------------
555 
556  kepler = false ;
557  for (int l=0; l<nzet; l++) {
558  int imax = mg->get_nr(l) - 1 ;
559  if (l == l_b) imax-- ; // The surface point is skipped
560  for (int i=0; i<imax; i++) {
561  if ( ent.val_grid_point(l, 0, j_b, i) < 0. ) {
562  kepler = true ;
563  cout << "ent < 0 for l, i : " << l << " " << i
564  << " ent = " << ent.val_grid_point(l, 0, j_b, i) << endl ;
565  }
566  }
567  }
568 
569  if ( kepler ) {
570  cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
571  omega_c /= fact_omega ; // Omega is decreased
572  cout << "New central rotation frequency : "
573  << omega_c/(2.*M_PI) * f_unit << " Hz" << endl ;
574  omega_trop_grand = true ;
575  }
576 
577  } // End of while ( kepler )
578 
579  if ( omega_trop_grand ) { // fact_omega is decreased for the
580  // next step
581  fact_omega = sqrt( fact_omega ) ;
582  cout << "**** New fact_omega : " << fact_omega << endl ;
583  }
584 
585  hydro_euler() ;
586 
587  beta = beta_new ;
588 
589  //---------------------------------------
590  // Calculate error of the GRV2 identity
591  //---------------------------------------
592 
593  err_grv2 = grv2() ;
594 
595 
596  //--------------------------------------
597  // Relaxations on some quantities....?
598  //
599  //--------------------------------------
600 
601  if (mer >= 10) {
602  logn = relax * logn + relax_prev * logn_prev ;
603 
604  psi = relax * psi_new + relax_prev * psi_prev ;
605 
606  beta = relax * beta + relax_prev * beta_prev ;
607 
608  }
609 
610  // Update of the metric quantities :
611 
612  update_metric() ;
613 
614  //-----------------------
615  // Informations display
616  // More to come later......
617  //-----------------------
618 
619  partial_display(cout) ;
620  fichfreq << " " << omega_c / (2*M_PI) * f_unit ;
621  fichevol << " " << ray_pole() / ray_eq() ;
622  fichevol << " " << ent_c ;
623 
624  //-----------------------------------------
625  // Convergence towards a given baryon mass
626  //-----------------------------------------
627 
628  if (mer > mer_mass) {
629 
630  double xx ;
631  if (mbar_wanted > 0.) {
632  xx = mass_b() / mbar_wanted - 1. ;
633  cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx
634  << endl ;
635  }
636  else{
637  xx = mass_g() / fabs(mbar_wanted) - 1. ;
638  cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx
639  << endl ;
640  }
641  double xprog = ( mer > 2*mer_mass) ? 1. :
642  double(mer-mer_mass)/double(mer_mass) ;
643  xx *= xprog ;
644  double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
645  double fact = pow(ax, aexp_mass) ;
646  cout << " xprog, xx, ax, fact : " << xprog << " " <<
647  xx << " " << ax << " " << fact << endl ;
648 
649  if ( change_ent ) {
650  ent_c *= fact ;
651  }
652  else {
653  if (mer%4 == 0) omega_c *= fact ;
654  }
655  }
656 
657 
658  //------------------------------------------------------------
659  // Relative change in enthalpy with respect to previous step
660  // ** Check: Is diffrel(ent, ent_prev) ok?
661  //------------------------------------------------------------
662 
663  Tbl diff_ent_tbl = diffrel( ent, ent_prev ) ;
664  diff_ent = diff_ent_tbl(0) ;
665  for (int l=1; l<nzet; l++) {
666  diff_ent += diff_ent_tbl(l) ;
667  }
668  diff_ent /= nzet ;
669 
670  fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
671  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
672  fichconv << " " << log10( fabs(max_triax) + 1.e-16 ) ;
673 
674  vit_triax = 0 ;
675  if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
676  vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
677  }
678 
679  fichconv << " " << vit_triax ;
680 
681  //------------------------------
682  // Recycling for the next step
683  //------------------------------
684 
685  ent_prev = ent ;
686  logn_prev = logn ;
687  psi_prev = psi ;
688  max_triax_prev = max_triax ;
689  beta_prev = beta ;
690 
691  fichconv << endl ;
692  fichfreq << endl ;
693  fichevol << endl ;
694  fichconv.flush() ;
695  fichfreq.flush() ;
696  fichevol.flush() ;
697 
698 
699  } // End of main loop
700 
701  //=========================================================================
702  // End of iteration
703  //=========================================================================
704 
705  fichconv.close() ;
706  fichfreq.close() ;
707  fichevol.close() ;
708 
709 }
710 }
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:676
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
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
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 .
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:151
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
Tbl par_frot
Parameters of the function .
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
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
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
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:97
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:630
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
Coord phi
coordinate centered on the grid
Definition: map.h:746
Coord sint
Definition: map.h:747
double ray_eq() const
Coordinate radius at , [r_unit].
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 homothetie(double lambda)=0
Sets a new radial scale.
double omega
Rotation angular velocity ([f_unit] )
Parameter storage.
Definition: param.h:125
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
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
void solve_psi(Scalar &psi_new)
Solution of the equations for the conformal factor for rotating stars in CFC.
Scalar psi4
Conformal factor .
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 .
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
virtual 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, Param *=0x0)
Computes an equilibrium configuration.
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
void fait_omega_field(double omeg_min, double omeg_max, double precis, int nitermax)
Computes (member omega_field ).
double omega_max
Maximum value of .
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition: cmp_math.C:325
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
int get_taille() const
Gives the total size (ie dim.taille)
Definition: tbl.h:417
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
double ray_pole() const
Coordinate radius at [r_unit].
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
double omega_min
Minimum value of .
Vector beta
Shift vector.
Definition: hot_star.h:133
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:616
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: hot_star.h:112