LORENE
strot_diff_equil.C
1 /*
2  * Function Star_rot_diff::equilibrium
3  *
4  * (see file star_rot_diff.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 "star_rot_diff.h"
33 #include "param.h"
34 #include "graphique.h"
35 #include "utilitaires.h"
36 #include "unites.h"
37 
38 #include "tenseur.h" // Unavoidable for now due to Poisson-solving functions
39 
40 
41 namespace Lorene{
42 void Star_rot_diff::equilibrium(double ent_c, double omega_c0, double fact_omega,
43  int nzadapt, const Tbl& ent_limit,
44  const Itbl& icontrol,
45  const Tbl& control, double mbar_wanted,
46  double aexp_mass, Tbl& diff, Param*){
47 
48  // Fundamental constants and units
49  // --------------------------------
50  using namespace Unites ;
51 
52  // For the display
53  // ---------------
54  char display_bold[]="x[1m" ; display_bold[0] = 27 ;
55  char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
56 
57 
58  // Grid parameters
59  // ----------------
60 
61  const Mg3d* mg = mp.get_mg() ;
62  int nz = mg->get_nzone() ; // total number of domains
63  int nzm1 = nz - 1 ;
64 
65  // The following is required to initialize mp_prev as a Map_et:
66  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
67 
68  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
69  assert(mg->get_type_t() == SYM) ;
70  int l_b = nzet - 1 ;
71  int i_b = mg->get_nr(l_b) - 1 ;
72  int j_b = mg->get_nt(l_b) - 1 ;
73  int k_b = 0 ;
74 
75  // Value of the enthalpy defining the surface of the star
76  double ent_b = ent_limit(nzet-1) ;
77 
78  // Parameters to control the iteration
79  // -----------------------------------
80 
81  int mer_max = icontrol(0) ;
82  int mer_rot = icontrol(1) ;
83  int mer_change_omega = icontrol(2) ;
84  int mer_fix_omega = icontrol(3) ;
85  int mer_mass = icontrol(4) ;
86  int mermax_poisson = icontrol(5) ;
87  int mer_triax = icontrol(6) ;
88  int delta_mer_kep = icontrol(7) ;
89 
90  // Protections:
91  if (mer_change_omega < mer_rot) {
92  cout << "Star_rot_diff::equilibrium: mer_change_omega < mer_rot !" << endl ;
93  cout << " mer_change_omega = " << mer_change_omega << endl ;
94  cout << " mer_rot = " << mer_rot << endl ;
95  abort() ;
96  }
97  if (mer_fix_omega < mer_change_omega) {
98  cout << "Star_rot_diff::equilibrium: mer_fix_omega < mer_change_omega !"
99  << endl ;
100  cout << " mer_fix_omega = " << mer_fix_omega << endl ;
101  cout << " mer_change_omega = " << mer_change_omega << endl ;
102  abort() ;
103  }
104 
105  // In order to converge to a given baryon mass, shall the central
106  // enthalpy be varied or Omega ?
107  bool change_ent = true ;
108  if (mer_mass < 0) {
109  change_ent = false ;
110  mer_mass = abs(mer_mass) ;
111  }
112 
113  double precis = control(0) ;
114  double omega_ini = control(1) ;
115  double relax = control(2) ;
116  double relax_prev = double(1) - relax ;
117  double relax_poisson = control(3) ;
118  double thres_adapt = control(4) ;
119  double ampli_triax = control(5) ;
120  double precis_adapt = control(6) ;
121 
122 
123  // Error indicators
124  // ----------------
125 
126  diff.set_etat_qcq() ;
127  double& diff_ent = diff.set(0) ;
128  double& diff_nuf = diff.set(1) ;
129  double& diff_nuq = diff.set(2) ;
130 // double& diff_dzeta = diff.set(3) ;
131 // double& diff_ggg = diff.set(4) ;
132  double& diff_shift_x = diff.set(5) ;
133  double& diff_shift_y = diff.set(6) ;
134  double& vit_triax = diff.set(7) ;
135 
136  // Parameters for the function Map_et::adapt
137  // -----------------------------------------
138 
139  Param par_adapt ;
140  int nitermax = 100 ;
141  int niter ;
142  int adapt_flag = 1 ; // 1 = performs the full computation,
143  // 0 = performs only the rescaling by
144  // the factor alpha_r
145  int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
146  // isosurfaces
147  double alpha_r ;
148  double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
149 
150  par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
151  // locate zeros by the secant method
152  par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment
153  // to the isosurfaces of ent is to be
154  // performed
155  par_adapt.add_int(nz_search, 2) ; // number of domains to search for
156  // the enthalpy isosurface
157  par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
158  // 0 = performs only the rescaling by
159  // the factor alpha_r
160  par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
161  // (theta_*, phi_*)
162  par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
163  // (theta_*, phi_*)
164 
165  par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
166  // the secant method
167 
168  par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in
169  // the determination of zeros by
170  // the secant method
171  par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping, 0 = contracting mapping
172 
173  par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
174  // distances will be multiplied
175 
176  par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
177  // to define the isosurfaces.
178 
179  // Parameters for the function Map_et::poisson for nuf
180  // ----------------------------------------------------
181 
182  double precis_poisson = 1.e-16 ;
183 
184  Cmp cssjm1_tggg(ssjm1_tggg) ;
185  Cmp cssjm1_khi(ssjm1_khi) ;
186  Tenseur cssjm1_wshift(mp, 1, CON, mp.get_bvect_cart() ) ;
187  cssjm1_wshift.set_etat_qcq() ;
188  for (int i=1; i<=3; i++) {
189  cssjm1_wshift.set(i-1) = ssjm1_wshift(i) ;
190  }
191 
192  Tenseur cshift(mp, 1, CON, mp.get_bvect_cart() ) ;
193  cshift.set_etat_qcq() ;
194  for (int i=1; i<=3; i++) {
195  cshift.set(i-1) = -beta(i) ;
196  }
197 
198  Tenseur cw_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
199  cw_shift.set_etat_qcq() ;
200  for (int i=1; i<=3; i++) {
201  cw_shift.set(i-1) = w_shift(i) ;
202  }
203 
204  Tenseur ckhi_shift(mp) ;
205  ckhi_shift.set_etat_qcq() ;
206  ckhi_shift.set() = khi_shift ;
207 
208  Param par_poisson_nuf ;
209  par_poisson_nuf.add_int(mermax_poisson, 0) ; // maximum number of iterations
210  par_poisson_nuf.add_double(relax_poisson, 0) ; // relaxation parameter
211  par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
212  par_poisson_nuf.add_int_mod(niter, 0) ; // number of iterations actually used
213  par_poisson_nuf.add_scalar_mod( ssjm1_nuf ) ;
214 
215  Param par_poisson_nuq ;
216  par_poisson_nuq.add_int(mermax_poisson, 0) ; // maximum number of iterations
217  par_poisson_nuq.add_double(relax_poisson, 0) ; // relaxation parameter
218  par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
219  par_poisson_nuq.add_int_mod(niter, 0) ; // number of iterations actually used
220  par_poisson_nuq.add_scalar_mod( ssjm1_nuq ) ;
221 
222  Param par_poisson_tggg ;
223  par_poisson_tggg.add_int(mermax_poisson, 0) ; // maximum number of iterations
224  par_poisson_tggg.add_double(relax_poisson, 0) ; // relaxation parameter
225  par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
226  par_poisson_tggg.add_int_mod(niter, 0) ; // number of iterations actually used
227  par_poisson_tggg.add_cmp_mod( cssjm1_tggg ) ;
228  double lambda_tggg ;
229  par_poisson_tggg.add_double_mod( lambda_tggg ) ;
230 
231  Param par_poisson_dzeta ;
232  double lbda_grv2 ;
233  par_poisson_dzeta.add_double_mod( lbda_grv2 ) ;
234 
235  // Parameters for the function Scalar::poisson_vect
236  // -------------------------------------------------
237 
238  Param par_poisson_vect ;
239 
240  par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
241  par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
242  par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
243  par_poisson_vect.add_cmp_mod( cssjm1_khi ) ;
244  par_poisson_vect.add_tenseur_mod( cssjm1_wshift ) ;
245  par_poisson_vect.add_int_mod(niter, 0) ;
246 
247 
248  // Initializations
249  // ---------------
250 
251  // Initial angular velocities
252  double omega_c = 0 ;
253 
254  double accrois_omega = (omega_c0 - omega_ini) /
255  double(mer_fix_omega - mer_change_omega) ;
256 
257 
258  update_metric() ; //update of the metric quantities
259 
260  equation_of_state() ; // update of the density, pressure,...etc
261 
262  hydro_euler() ; //update of the hydro quantities relative to the
263  // Eulerian observer
264 
265  // Quantities at the previous step :
266  Map_et mp_prev = mp_et ;
267  Scalar ent_prev = ent ;
268  Scalar logn_prev = logn ;
269  Scalar dzeta_prev = dzeta ;
270 
271  // Creation of uninitialized tensors:
272  Scalar source_nuf(mp) ; // source term in the equation for nuf
273  Scalar source_nuq(mp) ; // source term in the equation for nuq
274  Scalar source_dzf(mp) ; // matter source term in the eq. for dzeta
275  Scalar source_dzq(mp) ; // quadratic source term in the eq. for dzeta
276  Scalar source_tggg(mp) ; // source term in the eq. for tggg
277  Vector source_shift(mp, CON, mp.get_bvect_cart()) ; // source term for shift
278  Scalar mlngamma(mp) ; // centrifugal potential
279 
280  // Preparations for the Poisson equations:
281  // --------------------------------------
282  if (nuf.get_etat() == ETATZERO) {
283  nuf.set_etat_qcq() ;
284  nuf = 0 ;
285  }
286 
287  if (relativistic) {
288  if (nuq.get_etat() == ETATZERO) {
289  nuq.set_etat_qcq() ;
290  nuq = 0 ;
291  }
292 
293  if (tggg.get_etat() == ETATZERO) {
294  tggg.set_etat_qcq() ;
295  tggg = 0 ;
296  }
297 
298  if (dzeta.get_etat() == ETATZERO) {
299  dzeta.set_etat_qcq() ;
300  dzeta = 0 ;
301  }
302  }
303 
304  // Output files
305  // -------------
306 
307  ofstream fichconv("convergence.d") ; // Output file for diff_ent
308  fichconv << "# diff_ent GRV2 max_triax vit_triax" << endl ;
309 
310  ofstream fichfreq("frequency.d") ; // Output file for omega
311  fichfreq << "# f [Hz]" << endl ;
312 
313  ofstream fichevol("evolution.d") ; // Output file for various quantities
314  fichevol <<
315  "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
316  << endl ;
317 
318  diff_ent = 1 ;
319  double err_grv2 = 1 ;
320  double max_triax_prev = 0 ; // Triaxial amplitude at previous step
321 
322  //=========================================================================
323  // Start of iteration
324  //=========================================================================
325 
326  for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
327 
328  cout << "-----------------------------------------------" << endl ;
329  cout << "step: " << mer << endl ;
330  cout << "diff_ent = " << display_bold << diff_ent << display_normal
331  << endl ;
332  cout << "err_grv2 = " << err_grv2 << endl ;
333  fichconv << mer ;
334  fichfreq << mer ;
335  fichevol << mer ;
336 
337  // Switch on rotation
338  if (mer >= mer_rot) {
339 
340  if (mer < mer_change_omega) {
341  omega_c = omega_ini ;
342  }
343  else {
344  if (mer <= mer_fix_omega) {
345  omega_c = omega_ini + accrois_omega *
346  (mer - mer_change_omega) ;
347  }
348  }
349 
350  }
351 
352  //-----------------------------------------------
353  // Sources of the Poisson equations
354  //-----------------------------------------------
355 
356  // Source for nu
357  // -------------
358  Scalar bet = log(bbb) ;
359  bet.std_spectral_base() ;
360 
361  Vector d_logn = logn.derive_cov( mp.flat_met_spher() ) ;
362  Vector d_bet = bet.derive_cov( mp.flat_met_spher() ) ;
363 
364  if (relativistic) {
365  source_nuf = qpig * a_car *( ener_euler + s_euler ) ;
366 
367  source_nuq = ak_car - d_logn(1)*(d_logn(1)+d_bet(1))
368  - d_logn(2)*(d_logn(2)+d_bet(2))
369  - d_logn(3)*(d_logn(3)+d_bet(3)) ;
370  }
371  else {
372  source_nuf = qpig * nbar ;
373 
374  source_nuq = 0 ;
375  }
376  source_nuf.std_spectral_base() ;
377  source_nuq.std_spectral_base() ;
378 
379  // Source for dzeta
380  // ----------------
381  source_dzf = 2 * qpig * a_car * (press + (ener_euler+press) * uuu*uuu ) ;
382  source_dzf.std_spectral_base() ;
383 
384  source_dzq = 1.5 * ak_car
385  - d_logn(1)*d_logn(1) - d_logn(2)*d_logn(2) - d_logn(3)*d_logn(3) ;
386  source_dzq.std_spectral_base() ;
387 
388  // Source for tggg
389  // ---------------
390 
391  source_tggg = 4 * qpig * nn * a_car * bbb * press ;
392  source_tggg.std_spectral_base() ;
393 
394  source_tggg.mult_rsint() ;
395 
396 
397  // Source for shift
398  // ----------------
399 
400  // Matter term:
401  source_shift = (-4*qpig) * nn * a_car * (ener_euler + press)
402  * u_euler ;
403 
404  // Quadratic terms:
405  Vector vtmp = 6 * bet.derive_con( mp.flat_met_spher() )
406  - 2 * logn.derive_con( mp.flat_met_spher() ) ;
407  vtmp.change_triad(mp.get_bvect_cart()) ;
408 
409  Vector squad = nn * contract(tkij, 1, vtmp, 0) ;
410 
411  source_shift = source_shift + squad.up(0, mp.flat_met_cart() ) ;
412 
413  //----------------------------------------------
414  // Resolution of the Poisson equation for nuf
415  //----------------------------------------------
416 
417  source_nuf.poisson(par_poisson_nuf, nuf) ;
418 
419  cout << "Test of the Poisson equation for nuf :" << endl ;
420  Tbl err = source_nuf.test_poisson(nuf, cout, true) ;
421  diff_nuf = err(0, 0) ;
422 
423  //---------------------------------------
424  // Triaxial perturbation of nuf
425  //---------------------------------------
426 
427  if (mer == mer_triax) {
428 
429  if ( mg->get_np(0) == 1 ) {
430  cout <<
431  "Star_rot_diff::equilibrium: np must be stricly greater than 1"
432  << endl << " to set a triaxial perturbation !" << endl ;
433  abort() ;
434  }
435 
436  const Coord& phi = mp.phi ;
437  const Coord& sint = mp.sint ;
438  Scalar perturb(mp) ;
439  perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
440  nuf = nuf * perturb ;
441 
442  nuf.std_spectral_base() ; // set the bases for spectral expansions
443  // to be the standard ones for a
444  // scalar field
445 
446  }
447 
448  // Monitoring of the triaxial perturbation
449  // ---------------------------------------
450 
451  const Valeur& va_nuf = nuf.get_spectral_va() ;
452  va_nuf.coef() ; // Computes the spectral coefficients
453  double max_triax = 0 ;
454 
455  if ( mg->get_np(0) > 1 ) {
456 
457  for (int l=0; l<nz; l++) { // loop on the domains
458  for (int j=0; j<mg->get_nt(l); j++) {
459  for (int i=0; i<mg->get_nr(l); i++) {
460 
461  // Coefficient of cos(2 phi) :
462  double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ;
463 
464  // Coefficient of sin(2 phi) :
465  double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ;
466 
467  double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
468 
469  max_triax = ( xx > max_triax ) ? xx : max_triax ;
470  }
471  }
472  }
473 
474  }
475 
476  cout << "Triaxial part of nuf : " << max_triax << endl ;
477 
478  if (relativistic) {
479 
480  //----------------------------------------------
481  // Resolution of the Poisson equation for nuq
482  //----------------------------------------------
483 
484  source_nuq.poisson(par_poisson_nuq, nuq) ;
485 
486  cout << "Test of the Poisson equation for nuq :" << endl ;
487  err = source_nuq.test_poisson(nuq, cout, true) ;
488  diff_nuq = err(0, 0) ;
489 
490  //---------------------------------------------------------
491  // Resolution of the vector Poisson equation for the shift
492  //---------------------------------------------------------
493 
494  for (int i=1; i<=3; i++) {
495  if(source_shift(i).get_etat() != ETATZERO) {
496  if(source_shift(i).dz_nonzero()) {
497  assert( source_shift(i).get_dzpuis() == 4 ) ;
498  }
499  else{
500  (source_shift.set(i)).set_dzpuis(4) ;
501  }
502  }
503  }
504 
505  double lambda_shift = double(1) / double(3) ;
506 
507  if ( mg->get_np(0) == 1 ) {
508  lambda_shift = 0 ;
509  }
510 
511  Tenseur csource_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
512  csource_shift.set_etat_qcq() ;
513  if ( (source_shift(1).get_etat() == ETATZERO) &&
514  (source_shift(2).get_etat() == ETATZERO) &&
515  (source_shift(3).get_etat() == ETATZERO) ) {
516  csource_shift.set_etat_zero() ;
517  }
518  else {
519  csource_shift.set_etat_qcq() ;
520  for (int i=1; i<=3; i++) {
521  csource_shift.set(i-1) = source_shift(i) ;
522  }
523  csource_shift.set(2).set_etat_zero() ;
524  }
525 
526  csource_shift.poisson_vect(lambda_shift, par_poisson_vect,
527  cshift, cw_shift, ckhi_shift) ;
528  // beta = -source_shift.poisson(lambda_shift, par_poisson_vect, 2) ;
529 
530  for (int i=1; i<=3; i++) {
531  beta.set(i) = - cshift(i-1) ;
532  w_shift.set(i) = cw_shift(i-1) ;
533  }
534  khi_shift = ckhi_shift() ;
535 
536  cout << "Test of the Poisson equation for shift_x :" << endl ;
537  err = source_shift(1).test_poisson(-beta(1), cout, true) ;
538  diff_shift_x = err(0, 0) ;
539 
540  cout << "Test of the Poisson equation for shift_y :" << endl ;
541  err = source_shift(2).test_poisson(-beta(2), cout, true) ;
542  diff_shift_y = err(0, 0) ;
543 
544  // Computation of tnphi and nphi from the Cartesian components
545  // of the shift
546  // -----------------------------------------------------------
547 
548  fait_nphi() ;
549 
550  }
551 
552  //-----------------------------------------
553  // Determination of the fluid velociy U
554  //-----------------------------------------
555 
556  if (mer > mer_fix_omega + delta_mer_kep) {
557 
558  omega_c *= fact_omega ; // Increase of the angular velocity if
559  } // fact_omega != 1
560 
561 
562  bool omega_trop_grand = false ;
563  bool kepler = true ;
564 
565  while ( kepler ) {
566 
567  // Possible decrease of Omega to ensure a velocity < c
568 
569  bool superlum = true ;
570 
571  while ( superlum ) {
572 
573  // Computation of Omega(r,theta)
574 
575  if (omega_c == 0.) {
576  omega_field = 0 ;
577  }
578  else if (par_frot.get_taille() == 5) { // Uryu rotation profile
579  // First, compute F_e and F_max
580  Scalar brst2 = bbb ;
581  brst2.annule_domain(nzm1) ;
582  brst2.std_spectral_base() ;
583  brst2.mult_rsint() ; // Multiplication by r sin(theta)
584  brst2 = brst2*brst2 ;
585  Scalar nnn2 = nn * nn ;
586 
587  double omnp = omega_field.val_grid_point(l_b, k_b, j_b, i_b) - nphi.val_grid_point(l_b, k_b, j_b, i_b) ;
588  double F_e = brst2.val_grid_point(l_b, k_b, j_b, i_b) * omnp /
589  ( nnn2.val_grid_point(l_b, k_b, j_b, i_b) -
590  brst2.val_grid_point(l_b, k_b, j_b, i_b) * omnp * omnp ) ;
591 
592  if (F_e == double(0)) {} // Initial iteration
593  else {
594  int l, k, j, i ;
595  double om_max = 0 ;
596  double om_lkji, F_max ;
597  for (l=0; l<nzet+1; l++) {
598  for (k=0; k<mg->get_np(l); k++) {
599  for (j=0; j<mg->get_nt(l); j++) {
600  for (i=0; i<mg->get_nr(l); i++) {
601  om_lkji = omega_field.val_grid_point(l, k, j, i) ;
602  if (om_lkji > om_max) {
603  om_max = om_lkji ;
604  omnp = om_max - nphi.val_grid_point(l, k, j, i) ;
605  F_max = brst2.val_grid_point(l, k, j, i) * omnp /
606  ( nnn2.val_grid_point(l, k, j, i) -
607  brst2.val_grid_point(l, k, j, i) * omnp * omnp ) ;
608  }
609  }
610  }
611  }
612  }
613 
614  // Now A^2 and B^2, parameters for omega_frot, can be computed
615  double lambda1 = par_frot(1) ;
616  double lambda2 = par_frot(2) ;
617  int p = 1 ;
618  int q = 3 ;
619  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) )
620  / ((lambda1-1)*pow(F_e, p) - (lambda2-1)*pow(F_max, p)), 1./(p+q) ) ;
621  par_frot.set(4) = F_e*F_max/omega_c * pow((lambda2*pow(F_e, q)-lambda1*pow(F_max, q))
622  / (lambda2*(lambda1-1)*pow(F_e, p+q) - lambda1*(lambda2-1)*pow(F_max, p+q)), 1./p) ;
623 
624  cout << F_e << " " << F_max << endl ;
625  }
626 
627  par_frot.set(0) = omega_c ;
628 
629  cout << "Parameters of Omega(F) : " << par_frot << endl ;
630 
631  double omeg_min = 0 ;
632  double omeg_max = par_frot(1) * omega_c ;
633  double precis1 = 1.e-14 ;
634  int nitermax1 = 100 ;
635 
636  fait_omega_field(omeg_min, omeg_max, precis1, nitermax1) ;
637  }
638  else { // Usual rotation profile
639  par_frot.set(0) = omega_c ;
640  if (par_frot(2) != double(0)) { // fixed a = R_eq / R_0
641  par_frot.set(1) = ray_eq() / par_frot(2) ;
642  }
643  double omeg_min = 0 ;
644  double omeg_max = omega_c ;
645  double precis1 = 1.e-14 ;
646  int nitermax1 = 100 ;
647 
648  fait_omega_field(omeg_min, omeg_max, precis1, nitermax1) ;
649  }
650 
651  // New fluid velocity U :
652 
653  Scalar tmp = omega_field - nphi ;
654  tmp.annule_domain(nzm1) ;
655  tmp.std_spectral_base() ;
656 
657  tmp.mult_rsint() ; // Multiplication by r sin(theta)
658 
659  uuu = bbb / nn * tmp ;
660 
661  if (uuu.get_etat() == ETATQCQ) {
662  // Same basis as (Omega -N^phi) r sin(theta) :
663  (uuu.set_spectral_va()).set_base( tmp.get_spectral_va().get_base() ) ;
664  }
665 
666  // Is the new velocity larger than c in the equatorial plane ?
667 
668  superlum = false ;
669 
670  for (int l=0; l<nzet; l++) {
671  for (int i=0; i<mg->get_nr(l); i++) {
672 
673  double u1 = uuu.val_grid_point(l, 0, j_b, i) ;
674  if (u1 >= 1.) { // superluminal velocity
675  superlum = true ;
676  cout << "U > c for l, i : " << l << " " << i
677  << " U = " << u1 << endl ;
678  }
679  }
680  }
681  if ( superlum ) {
682  cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ;
683  omega_c /= fact_omega ; // Decrease of Omega_c
684  cout << "New central rotation frequency : "
685  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
686  omega_trop_grand = true ;
687  }
688  } // end of while ( superlum )
689 
690 
691  // New computation of U (which this time is not superluminal)
692  // as well as of gam_euler, ener_euler, etc...
693  // -----------------------------------
694 
695  hydro_euler() ;
696 
697 
698  //------------------------------------------------------
699  // First integral of motion
700  //------------------------------------------------------
701 
702  // Centrifugal potential :
703  if (relativistic) {
704  mlngamma = - log( gam_euler ) ;
705  }
706  else {
707  mlngamma = - 0.5 * uuu*uuu ;
708  }
709 
710  // Equatorial values of various potentials :
711  double nuf_b = nuf.val_grid_point(l_b, k_b, j_b, i_b) ;
712  double nuq_b = nuq.val_grid_point(l_b, k_b, j_b, i_b) ;
713  double mlngamma_b = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ;
714  double primf_b = prim_field.val_grid_point(l_b, k_b, j_b, i_b) ;
715 
716 
717  // Central values of various potentials :
718  double nuf_c = nuf.val_grid_point(0,0,0,0) ;
719  double nuq_c = nuq.val_grid_point(0,0,0,0) ;
720  double mlngamma_c = 0 ;
721  double primf_c = prim_field.val_grid_point(0,0,0,0) ;
722 
723  // Scale factor to ensure that the enthalpy is equal to ent_b at
724  // the equator
725  double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
726  + nuq_c - nuq_b + primf_c - primf_b)
727  / ( nuf_b - nuf_c ) ;
728  alpha_r = sqrt(alpha_r2) ;
729  cout << "alpha_r = " << alpha_r << endl ;
730 
731  // Readjustment of nu :
732  // -------------------
733 
734  logn = alpha_r2 * nuf + nuq ;
735  double nu_c = logn.val_grid_point(0,0,0,0) ;
736 
737  // First integral --> enthalpy in all space
738  //-----------------
739 
740  ent = (ent_c + nu_c + mlngamma_c) - logn - mlngamma - prim_field ;
741 
742  // Test: is the enthalpy negative somewhere in the equatorial plane
743  // inside the star ? If yes, this means that the Keplerian velocity
744  // has been overstep.
745 
746  kepler = false ;
747  for (int l=0; l<nzet; l++) {
748  int imax = mg->get_nr(l) - 1 ;
749  if (l == l_b) imax-- ; // The surface point is skipped
750  for (int i=0; i<imax; i++) {
751  if ( ent.val_grid_point(l, 0, j_b, i) < 0. ) {
752  kepler = true ;
753  cout << "ent < 0 for l, i : " << l << " " << i
754  << " ent = " << ent.val_grid_point(l, 0, j_b, i) << endl ;
755  }
756  }
757  }
758 
759  if ( kepler ) {
760  cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
761  omega_c /= fact_omega ; // Omega is decreased
762  cout << "New central rotation frequency : "
763  << omega_c/(2.*M_PI) * f_unit << " Hz" << endl ;
764  omega_trop_grand = true ;
765  }
766 
767  } // End of while ( kepler )
768 
769  if ( omega_trop_grand ) { // fact_omega is decreased for the
770  // next step
771  fact_omega = sqrt( fact_omega ) ;
772  cout << "**** New fact_omega : " << fact_omega << endl ;
773  }
774 
775  //----------------------------------------------------
776  // Adaptation of the mapping to the new enthalpy field
777  //----------------------------------------------------
778 
779  // Shall the adaptation be performed (cusp) ?
780  // ------------------------------------------
781 
782  double dent_eq = ent.dsdr().val_grid_point(l_b, k_b, j_b, i_b) ;
783  double dent_pole = ent.dsdr().val_grid_point(l_b, k_b, 0, i_b) ;
784  double rap_dent = fabs( dent_eq / dent_pole ) ;
785  cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
786 
787  if ( rap_dent < thres_adapt ) {
788  adapt_flag = 0 ; // No adaptation of the mapping
789  cout << "******* FROZEN MAPPING *********" << endl ;
790  }
791  else{
792  adapt_flag = 1 ; // The adaptation of the mapping is to be
793  // performed
794  }
795 
796  mp_prev = mp_et ;
797 
798  Cmp cent(ent) ;
799  mp.adapt(cent, par_adapt) ;
800 
801  //----------------------------------------------------
802  // Computation of the enthalpy at the new grid points
803  //----------------------------------------------------
804 
805  mp_prev.homothetie(alpha_r) ;
806 
807  mp.reevaluate(&mp_prev, nzet+1, cent) ;
808  ent = cent ;
809 
810  //----------------------------------------------------
811  // Equation of state
812  //----------------------------------------------------
813 
814  equation_of_state() ; // computes new values for nbar (n), ener (e)
815  // and press (p) from the new ent (H)
816 
817  //---------------------------------------------------------
818  // Matter source terms in the gravitational field equations
819  //---------------------------------------------------------
820 
821  //## Computation of tnphi and nphi from the Cartesian components
822  // of the shift for the test in hydro_euler():
823 
824  fait_nphi() ;
825 
826  hydro_euler() ; // computes new values for ener_euler (E),
827  // s_euler (S) and u_euler (U^i)
828 
829  if (relativistic) {
830 
831  //-------------------------------------------------------
832  // 2-D Poisson equation for tggg
833  //-------------------------------------------------------
834 
835  Cmp csource_tggg(source_tggg) ;
836  Cmp ctggg(tggg) ;
837  mp.poisson2d(csource_tggg, mp.cmp_zero(), par_poisson_tggg,
838  ctggg) ;
839  tggg = ctggg ;
840 
841 
842  //-------------------------------------------------------
843  // 2-D Poisson equation for dzeta
844  //-------------------------------------------------------
845 
846  Cmp csource_dzf(source_dzf) ;
847  Cmp csource_dzq(source_dzq) ;
848  Cmp cdzeta(dzeta) ;
849  mp.poisson2d(csource_dzf, csource_dzq, par_poisson_dzeta,
850  cdzeta) ;
851  dzeta = cdzeta ;
852 
853  err_grv2 = lbda_grv2 - 1;
854  cout << "GRV2: " << err_grv2 << endl ;
855 
856  }
857  else {
858  err_grv2 = grv2() ;
859  }
860 
861 
862  //---------------------------------------
863  // Computation of the metric coefficients (except for N^phi)
864  //---------------------------------------
865 
866  // Relaxations on nu and dzeta :
867 
868  if (mer >= 10) {
869  logn = relax * logn + relax_prev * logn_prev ;
870 
871  dzeta = relax * dzeta + relax_prev * dzeta_prev ;
872  }
873 
874  // Update of the metric coefficients N, A, B and computation of K_ij :
875 
876  update_metric() ;
877 
878  //-----------------------
879  // Informations display
880  //-----------------------
881 
882  partial_display(cout) ;
883  fichfreq << " " << omega_c / (2*M_PI) * f_unit ;
884  fichevol << " " << rap_dent ;
885  fichevol << " " << ray_pole() / ray_eq() ;
886  fichevol << " " << ent_c ;
887 
888  //-----------------------------------------
889  // Convergence towards a given baryon mass
890  //-----------------------------------------
891 
892  if (mer > mer_mass) {
893 
894  double xx ;
895  if (mbar_wanted > 0.) {
896  xx = mass_b() / mbar_wanted - 1. ;
897  cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx
898  << endl ;
899  }
900  else{
901  xx = mass_g() / fabs(mbar_wanted) - 1. ;
902  cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx
903  << endl ;
904  }
905  double xprog = ( mer > 2*mer_mass) ? 1. :
906  double(mer-mer_mass)/double(mer_mass) ;
907  xx *= xprog ;
908  double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
909  double fact = pow(ax, aexp_mass) ;
910  cout << " xprog, xx, ax, fact : " << xprog << " " <<
911  xx << " " << ax << " " << fact << endl ;
912 
913  if ( change_ent ) {
914  ent_c *= fact ;
915  }
916  else {
917  if (mer%4 == 0) omega_c *= fact ;
918  }
919  }
920 
921 
922  //------------------------------------------------------------
923  // Relative change in enthalpy with respect to previous step
924  //------------------------------------------------------------
925 
926  Tbl diff_ent_tbl = diffrel( ent, ent_prev ) ;
927  diff_ent = diff_ent_tbl(0) ;
928  for (int l=1; l<nzet; l++) {
929  diff_ent += diff_ent_tbl(l) ;
930  }
931  diff_ent /= nzet ;
932 
933  fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
934  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
935  fichconv << " " << log10( fabs(max_triax) + 1.e-16 ) ;
936 
937  vit_triax = 0 ;
938  if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
939  vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
940  }
941 
942  fichconv << " " << vit_triax ;
943 
944  //------------------------------
945  // Recycling for the next step
946  //------------------------------
947 
948  ent_prev = ent ;
949  logn_prev = logn ;
950  dzeta_prev = dzeta ;
951  max_triax_prev = max_triax ;
952 
953  fichconv << endl ;
954  fichfreq << endl ;
955  fichevol << endl ;
956  fichconv.flush() ;
957  fichfreq.flush() ;
958  fichevol.flush() ;
959 
960  } // End of main loop
961 
962  //=========================================================================
963  // End of iteration
964  //=========================================================================
965 
966  ssjm1_tggg = cssjm1_tggg ;
967  ssjm1_khi = cssjm1_khi ;
968  for (int i=1; i<=3; i++) {
969  ssjm1_wshift.set(i) = cssjm1_wshift(i-1) ;
970  }
971 
972  fichconv.close() ;
973  fichfreq.close() ;
974  fichevol.close() ;
975 
976 }
977 }
virtual double mass_g() const
Gravitational mass.
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:676
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
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition: param.C:1145
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: star_rot.h:237
Scalar dzeta
Metric potential .
Definition: star_rot.h:146
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:116
Radial mapping of rather general form.
Definition: map.h:2859
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition: param.C:249
Scalar prim_field
Field .
Definition: star_rot_diff.h:79
Map & mp
Mapping associated with the star.
Definition: star.h:180
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: star_rot.h:138
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
Definition: map.h:834
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
Tbl par_frot
Parameters of the function .
Definition: star_rot_diff.h:70
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: star_rot.h:172
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: star_rot.h:106
Scalar bbb
Metric factor B.
Definition: star_rot.h:119
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
Tensor up(int ind, const Metric &gam) const
Computes a new tensor by raising an index of *this.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition: scalar_pde.C:139
virtual double grv2() const
Error on the virial identity GRV2.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
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
void update_metric()
Computes metric coefficients from known potentials.
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
Scalar nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: star_rot.h:143
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot.h:113
Scalar nbar
Baryon density in the fluid frame.
Definition: star.h:192
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:566
Scalar ent
Log-enthalpy.
Definition: star.h:190
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition: star.h:201
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
Vector beta
Shift vector.
Definition: star.h:228
Tbl test_poisson(const Scalar &uu, ostream &ostr, bool detail=false) const
Checks if a Poisson equation with *this as a source has been correctly solved.
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Definition: scalar_deriv.C:401
Tensor field of valence 1.
Definition: vector.h:188
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:97
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:337
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Definition: star_rot.h:220
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
Definition: param.C:456
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:649
Coord phi
coordinate centered on the grid
Definition: map.h:747
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
Coord sint
Definition: map.h:748
Scalar nphi
Metric coefficient .
Definition: star_rot.h:125
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.
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: star.h:204
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
virtual void poisson2d(const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const =0
Computes the solution of a 2-D Poisson equation.
Scalar press
Fluid pressure.
Definition: star.h:194
virtual double mass_b() const
Baryon mass.
Scalar ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition: star_rot.h:210
Parameter storage.
Definition: param.h:125
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition: param.C:525
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition: map_et.C:928
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: star.h:207
Sym_tensor tkij
Tensor related to the extrinsic curvature tensor by .
Definition: star_rot.h:179
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 .
Scalar logn
Logarithm of the lapse N .
Definition: star.h:222
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
double ray_pole() const
Coordinate radius at [r_unit].
Definition: star_global.C:281
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
Multi-domain grid.
Definition: grilles.h:276
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition: star_rot.C:613
void add_scalar_mod(Scalar &ti, int position=0)
Adds the address of a new modifiable Scalar to the list.
Definition: param.C:1420
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:818
Scalar nn
Lapse function N .
Definition: star.h:225
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition: cmp_math.C:325
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: star_rot.h:228
int get_taille() const
Gives the total size (ie dim.taille)
Definition: tbl.h:417
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:133
virtual void reevaluate(const Map *mp_prev, int nzet, Cmp &uu) const =0
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:318
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:112
Scalar ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition: star_rot.h:204
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
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:616
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition: param.C:1007
Scalar tggg
Metric potential .
Definition: star_rot.h:149
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: star.C:465
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Definition: star_rot.C:779
Scalar ak_car
Scalar .
Definition: star_rot.h:198
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:389
Vector w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: star_rot.h:162
void fait_omega_field(double omeg_min, double omeg_max, double precis, int nitermax)
Computes (member omega_field ).
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition: star.h:198
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:327
Scalar omega_field
Field .
Definition: star_rot_diff.h:73
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition: param.C:388
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:613