LORENE
star_rot_equil.C
1 /*
2  * Method Star_rot::equilibrium
3  *
4  * (see file star_h.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 2010 Eric Gourgoulhon
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 as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 
31 
32 /*
33  * $Id: star_rot_equil.C,v 1.11 2025/03/06 09:42:41 j_novak Exp $
34  * $Log: star_rot_equil.C,v $
35  * Revision 1.11 2025/03/06 09:42:41 j_novak
36  * Definition of Map_XX::integrale(Scalar&).
37  *
38  * Revision 1.10 2025/03/04 13:16:50 j_novak
39  * New complete versions of Map_af::poisson() and Map_et::poisson() for Scalar, not using Cmp.
40  *
41  * Revision 1.9 2024/10/25 15:33:43 j_novak
42  * Added an exception handling if Kepler frequency is reached and we are not searching it (experimental).
43  *
44  * Revision 1.8 2022/02/14 10:50:24 j_novak
45  * Increased the number of digits in outputting files.
46  *
47  * Revision 1.7 2016/12/05 16:18:15 j_novak
48  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
49  *
50  * Revision 1.6 2014/10/13 08:53:39 j_novak
51  * Lorene classes and functions now belong to the namespace Lorene.
52  *
53  * Revision 1.5 2014/10/06 15:13:17 j_novak
54  * Modified #include directives to use c++ syntax.
55  *
56  * Revision 1.4 2010/01/26 16:49:53 e_gourgoulhon
57  * Reformated some outputs to the screen.
58  *
59  * Revision 1.3 2010/01/26 10:49:43 e_gourgoulhon
60  * First working version.
61  *
62  * Revision 1.2 2010/01/25 22:32:38 e_gourgoulhon
63  * Start of real implementation
64  *
65  * Revision 1.1 2010/01/25 18:15:52 e_gourgoulhon
66  * First version.
67  *
68  *
69  * $Header: /cvsroot/Lorene/C++/Source/Star/star_rot_equil.C,v 1.11 2025/03/06 09:42:41 j_novak Exp $
70  *
71  */
72 
73 // Headers Lorene
74 #include "star_rot.h"
75 #include "param.h"
76 #include "tenseur.h"
77 
78 #include "utilitaires.h"
79 #include "unites.h"
80 
81 namespace Lorene {
82 void Star_rot::equilibrium(double ent_c, double omega0, double fact_omega,
83  int nzadapt, const Tbl& ent_limit, const Itbl& icontrol,
84  const Tbl& control, double mbar_wanted,
85  double aexp_mass, Tbl& diff, Param*) {
86 
87  // Fundamental constants and units
88  // -------------------------------
89 
90  using namespace Unites ;
91 
92  // For the display
93  // ---------------
94  char display_bold[]="x[1m" ; display_bold[0] = 27 ;
95  char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
96 
97  // Grid parameters
98  // ---------------
99 
100  const Mg3d* mg = mp.get_mg() ;
101  int nz = mg->get_nzone() ; // total number of domains
102  int nzm1 = nz - 1 ;
103 
104  // The following is required to initialize mp_prev as a Map_et:
105  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
106 
107  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
108  assert(mg->get_type_t() == SYM) ;
109  int l_b = nzet - 1 ;
110  int i_b = mg->get_nr(l_b) - 1 ;
111  int j_b = mg->get_nt(l_b) - 1 ;
112  int k_b = 0 ;
113 
114  // Value of the enthalpy defining the surface of the star
115  double ent_b = ent_limit(nzet-1) ;
116 
117  // Parameters to control the iteration
118  // -----------------------------------
119 
120  int mer_max = icontrol(0) ;
121  int mer_rot = icontrol(1) ;
122  int mer_change_omega = icontrol(2) ;
123  int mer_fix_omega = icontrol(3) ;
124  int mer_mass = icontrol(4) ;
125  int mermax_poisson = icontrol(5) ;
126  int mer_triax = icontrol(6) ;
127  int delta_mer_kep = icontrol(7) ;
128 
129  // Protections:
130  if (mer_change_omega < mer_rot) {
131  cout << "Star_rot::equilibrium: mer_change_omega < mer_rot !" << endl ;
132  cout << " mer_change_omega = " << mer_change_omega << endl ;
133  cout << " mer_rot = " << mer_rot << endl ;
134  abort() ;
135  }
136  if (mer_fix_omega < mer_change_omega) {
137  cout << "Star_rot::equilibrium: mer_fix_omega < mer_change_omega !"
138  << endl ;
139  cout << " mer_fix_omega = " << mer_fix_omega << endl ;
140  cout << " mer_change_omega = " << mer_change_omega << endl ;
141  abort() ;
142  }
143 
144  // In order to converge to a given baryon mass, shall the central
145  // enthalpy be varied or Omega ?
146  bool change_ent = true ;
147  if (mer_mass < 0) {
148  change_ent = false ;
149  mer_mass = abs(mer_mass) ;
150  }
151 
152  double precis = control(0) ;
153  double omega_ini = control(1) ;
154  double relax = control(2) ;
155  double relax_prev = double(1) - relax ;
156  double relax_poisson = control(3) ;
157  double thres_adapt = control(4) ;
158  double ampli_triax = control(5) ;
159  double precis_adapt = control(6) ;
160 
161 
162  // Error indicators
163  // ----------------
164 
165  diff.set_etat_qcq() ;
166  double& diff_ent = diff.set(0) ;
167  double& diff_nuf = diff.set(1) ;
168  double& diff_nuq = diff.set(2) ;
169 // double& diff_dzeta = diff.set(3) ;
170 // double& diff_ggg = diff.set(4) ;
171  double& diff_shift_x = diff.set(5) ;
172  double& diff_shift_y = diff.set(6) ;
173  double& vit_triax = diff.set(7) ;
174 
175  // Parameters for the function Map_et::adapt
176  // -----------------------------------------
177 
178  Param par_adapt ;
179  int nitermax = 100 ;
180  int niter ;
181  int adapt_flag = 1 ; // 1 = performs the full computation,
182  // 0 = performs only the rescaling by
183  // the factor alpha_r
184  int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
185  // isosurfaces
186  double alpha_r ;
187  double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
188 
189  par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
190  // locate zeros by the secant method
191  par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment
192  // to the isosurfaces of ent is to be
193  // performed
194  par_adapt.add_int(nz_search, 2) ; // number of domains to search for
195  // the enthalpy isosurface
196  par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
197  // 0 = performs only the rescaling by
198  // the factor alpha_r
199  par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
200  // (theta_*, phi_*)
201  par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
202  // (theta_*, phi_*)
203 
204  par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
205  // the secant method
206 
207  par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in
208  // the determination of zeros by
209  // the secant method
210  par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping, 0 = contracting mapping
211 
212  par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
213  // distances will be multiplied
214 
215  par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
216  // to define the isosurfaces.
217 
218  // Parameters for the function Map_et::poisson for nuf
219  // ----------------------------------------------------
220 
221  double precis_poisson = 1.e-16 ;
222 
223  Cmp cssjm1_tggg(ssjm1_tggg) ;
224  Cmp cssjm1_khi(ssjm1_khi) ;
225  Tenseur cssjm1_wshift(mp, 1, CON, mp.get_bvect_cart() ) ;
226  cssjm1_wshift.set_etat_qcq() ;
227  for (int i=1; i<=3; i++) {
228  cssjm1_wshift.set(i-1) = ssjm1_wshift(i) ;
229  }
230 
231  Tenseur cshift(mp, 1, CON, mp.get_bvect_cart() ) ;
232  cshift.set_etat_qcq() ;
233  for (int i=1; i<=3; i++) {
234  cshift.set(i-1) = -beta(i) ;
235  }
236 
237  Tenseur cw_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
238  cw_shift.set_etat_qcq() ;
239  for (int i=1; i<=3; i++) {
240  cw_shift.set(i-1) = w_shift(i) ;
241  }
242 
243  Tenseur ckhi_shift(mp) ;
244  ckhi_shift.set_etat_qcq() ;
245  ckhi_shift.set() = khi_shift ;
246 
247  Param par_poisson_nuf ;
248  par_poisson_nuf.add_int(mermax_poisson, 0) ; // maximum number of iterations
249  par_poisson_nuf.add_double(relax_poisson, 0) ; // relaxation parameter
250  par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
251  par_poisson_nuf.add_int_mod(niter, 0) ; // number of iterations actually used
252  par_poisson_nuf.add_scalar_mod( ssjm1_nuf ) ;
253 
254  Param par_poisson_nuq ;
255  par_poisson_nuq.add_int(mermax_poisson, 0) ; // maximum number of iterations
256  par_poisson_nuq.add_double(relax_poisson, 0) ; // relaxation parameter
257  par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
258  par_poisson_nuq.add_int_mod(niter, 0) ; // number of iterations actually used
259  par_poisson_nuq.add_scalar_mod( ssjm1_nuq ) ;
260 
261  Param par_poisson_tggg ;
262  par_poisson_tggg.add_int(mermax_poisson, 0) ; // maximum number of iterations
263  par_poisson_tggg.add_double(relax_poisson, 0) ; // relaxation parameter
264  par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
265  par_poisson_tggg.add_int_mod(niter, 0) ; // number of iterations actually used
266  par_poisson_tggg.add_cmp_mod( cssjm1_tggg ) ;
267  double lambda_tggg ;
268  par_poisson_tggg.add_double_mod( lambda_tggg ) ;
269 
270  Param par_poisson_dzeta ;
271  double lbda_grv2 ;
272  par_poisson_dzeta.add_double_mod( lbda_grv2 ) ;
273 
274  // Parameters for the function Scalar::poisson_vect
275  // -------------------------------------------------
276 
277  Param par_poisson_vect ;
278 
279  par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
280  par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
281  par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
282  par_poisson_vect.add_cmp_mod( cssjm1_khi ) ;
283  par_poisson_vect.add_tenseur_mod( cssjm1_wshift ) ;
284  par_poisson_vect.add_int_mod(niter, 0) ;
285 
286 
287  // Initializations
288  // ---------------
289 
290  // Initial angular velocity
291  omega = 0 ;
292 
293  double accrois_omega = (omega0 - omega_ini) /
294  double(mer_fix_omega - mer_change_omega) ;
295 
296 
297  update_metric() ; // update of the metric coefficients
298 
299  equation_of_state() ; // update of the density, pressure, etc...
300 
301  hydro_euler() ; // update of the hydro quantities relative to the
302  // Eulerian observer
303 
304  // Quantities at the previous step :
305  Map_et mp_prev = mp_et ;
306  Scalar ent_prev = ent ;
307  Scalar logn_prev = logn ;
308  Scalar dzeta_prev = dzeta ;
309 
310  // Creation of uninitialized tensors:
311  Scalar source_nuf(mp) ; // source term in the equation for nuf
312  Scalar source_nuq(mp) ; // source term in the equation for nuq
313  Scalar source_dzf(mp) ; // matter source term in the eq. for dzeta
314  Scalar source_dzq(mp) ; // quadratic source term in the eq. for dzeta
315  Scalar source_tggg(mp) ; // source term in the eq. for tggg
316  Vector source_shift(mp, CON, mp.get_bvect_cart()) ;
317  // source term for shift
318  Scalar mlngamma(mp) ; // centrifugal potential
319 
320 
321  ofstream fichconv("convergence.d") ; // Output file for diff_ent
322  fichconv << "# diff_ent GRV2 max_triax vit_triax" << endl ;
323 
324  ofstream fichfreq("frequency.d") ; // Output file for omega
325  fichfreq << "# f [Hz]" << endl ;
326 
327  ofstream fichevol("evolution.d") ; // Output file for various quantities
328  fichevol <<
329  "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
330  << endl ;
331 
332  diff_ent = 1 ;
333  double err_grv2 = 1 ;
334  double max_triax_prev = 0 ; // Triaxial amplitude at previous step
335 
336  //=========================================================================
337  // Start of iteration
338  //=========================================================================
339 
340  for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
341 
342  cout << "-----------------------------------------------" << endl ;
343  cout << "step: " << mer << endl ;
344  cout << "diff_ent = " << display_bold << diff_ent << display_normal
345  << endl ;
346  cout << "err_grv2 = " << err_grv2 << endl ;
347  fichconv << mer ;
348  fichfreq << mer ;
349  fichevol << mer ;
350 
351  if (mer >= mer_rot) {
352 
353  if (mer < mer_change_omega) {
354  omega = omega_ini ;
355  }
356  else {
357  if (mer <= mer_fix_omega) {
358  omega = omega_ini + accrois_omega *
359  (mer - mer_change_omega) ;
360  }
361  }
362 
363  }
364 
365  //-----------------------------------------------
366  // Sources of the Poisson equations
367  //-----------------------------------------------
368 
369  // Source for nu
370  // -------------
371  Scalar bet = log(bbb) ;
372  bet.std_spectral_base() ;
373 
374  Vector d_logn = logn.derive_cov( mp.flat_met_spher() ) ;
375  Vector d_bet = bet.derive_cov( mp.flat_met_spher() ) ;
376 
377  if (relativistic) {
378  source_nuf = qpig * a_car *( ener_euler + s_euler ) ;
379 
380  source_nuq = ak_car - d_logn(1)*(d_logn(1)+d_bet(1))
381  - d_logn(2)*(d_logn(2)+d_bet(2))
382  - d_logn(3)*(d_logn(3)+d_bet(3)) ;
383  }
384  else {
385  source_nuf = qpig * nbar ;
386 
387  source_nuq = 0 ;
388  }
389  source_nuf.std_spectral_base() ;
390  source_nuq.std_spectral_base() ;
391 
392  // Source for dzeta
393  // ----------------
394  source_dzf = 2 * qpig * a_car * (press + (ener_euler+press) * uuu*uuu ) ;
395  source_dzf.std_spectral_base() ;
396 
397  source_dzq = 1.5 * ak_car
398  - d_logn(1)*d_logn(1) - d_logn(2)*d_logn(2) - d_logn(3)*d_logn(3) ;
399  source_dzq.std_spectral_base() ;
400 
401  // Source for tggg
402  // ---------------
403 
404  source_tggg = 4 * qpig * nn * a_car * bbb * press ;
405  source_tggg.std_spectral_base() ;
406 
407  source_tggg.mult_rsint() ;
408 
409 
410  // Source for shift
411  // ----------------
412 
413  // Matter term:
414  source_shift = (-4*qpig) * nn * a_car * (ener_euler + press)
415  * u_euler ;
416 
417  // Quadratic terms:
418  Vector vtmp = 6 * bet.derive_con( mp.flat_met_spher() )
419  - 2 * logn.derive_con( mp.flat_met_spher() ) ;
420  vtmp.change_triad(mp.get_bvect_cart()) ;
421 
422  Vector squad = nn * contract(tkij, 1, vtmp, 0) ;
423 
424  source_shift = source_shift + squad.up(0, mp.flat_met_cart() ) ;
425 
426  //----------------------------------------------
427  // Resolution of the Poisson equation for nuf
428  //----------------------------------------------
429 
430  source_nuf.poisson(par_poisson_nuf, nuf) ;
431 
432  cout << "Test of the Poisson equation for nuf :" << endl ;
433  Tbl err = source_nuf.test_poisson(nuf, cout, true) ;
434  diff_nuf = err(0, 0) ;
435 
436  //---------------------------------------
437  // Triaxial perturbation of nuf
438  //---------------------------------------
439 
440  if (mer == mer_triax) {
441 
442  if ( mg->get_np(0) == 1 ) {
443  cout <<
444  "Star_rot::equilibrium: np must be stricly greater than 1"
445  << endl << " to set a triaxial perturbation !" << endl ;
446  abort() ;
447  }
448 
449  const Coord& phi = mp.phi ;
450  const Coord& sint = mp.sint ;
451  Scalar perturb(mp) ;
452  perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
453  nuf = nuf * perturb ;
454 
455  nuf.std_spectral_base() ; // set the bases for spectral expansions
456  // to be the standard ones for a
457  // scalar field
458 
459  }
460 
461  // Monitoring of the triaxial perturbation
462  // ---------------------------------------
463 
464  const Valeur& va_nuf = nuf.get_spectral_va() ;
465  va_nuf.coef() ; // Computes the spectral coefficients
466  double max_triax = 0 ;
467 
468  if ( mg->get_np(0) > 1 ) {
469 
470  for (int l=0; l<nz; l++) { // loop on the domains
471  for (int j=0; j<mg->get_nt(l); j++) {
472  for (int i=0; i<mg->get_nr(l); i++) {
473 
474  // Coefficient of cos(2 phi) :
475  double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ;
476 
477  // Coefficient of sin(2 phi) :
478  double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ;
479 
480  double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
481 
482  max_triax = ( xx > max_triax ) ? xx : max_triax ;
483  }
484  }
485  }
486 
487  }
488 
489  cout << "Triaxial part of nuf : " << max_triax << endl ;
490 
491  if (relativistic) {
492 
493  //----------------------------------------------
494  // Resolution of the Poisson equation for nuq
495  //----------------------------------------------
496 
497  source_nuq.poisson(par_poisson_nuq, nuq) ;
498 
499  cout << "Test of the Poisson equation for nuq :" << endl ;
500  err = source_nuq.test_poisson(nuq, cout, true) ;
501  diff_nuq = err(0, 0) ;
502 
503  //---------------------------------------------------------
504  // Resolution of the vector Poisson equation for the shift
505  //---------------------------------------------------------
506 
507 
508  for (int i=1; i<=3; i++) {
509  if(source_shift(i).get_etat() != ETATZERO) {
510  if(source_shift(i).dz_nonzero()) {
511  assert( source_shift(i).get_dzpuis() == 4 ) ;
512  }
513  else{
514  (source_shift.set(i)).set_dzpuis(4) ;
515  }
516  }
517  }
518 
519  double lambda_shift = double(1) / double(3) ;
520 
521  if ( mg->get_np(0) == 1 ) {
522  lambda_shift = 0 ;
523  }
524 
525  Tenseur csource_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
526  csource_shift.set_etat_qcq() ;
527  for (int i=1; i<=3; i++) {
528  csource_shift.set(i-1) = source_shift(i) ;
529  }
530  csource_shift.set(2).set_etat_zero() ; //## bizarre...
531 
532  csource_shift.poisson_vect(lambda_shift, par_poisson_vect,
533  cshift, cw_shift, ckhi_shift) ;
534 
535  for (int i=1; i<=3; i++) {
536  beta.set(i) = - cshift(i-1) ;
537  beta.set(i).set_dzpuis(0) ; //## bizarre...
538  w_shift.set(i) = cw_shift(i-1) ;
539  }
540  khi_shift = ckhi_shift() ;
541 
542  cout << "Test of the Poisson equation for shift_x :" << endl ;
543  err = source_shift(1).test_poisson(-beta(1), cout, true) ;
544  diff_shift_x = err(0, 0) ;
545 
546  cout << "Test of the Poisson equation for shift_y :" << endl ;
547  err = source_shift(2).test_poisson(-beta(2), cout, true) ;
548  diff_shift_y = err(0, 0) ;
549 
550  // Computation of tnphi and nphi from the Cartesian components
551  // of the shift
552  // -----------------------------------------------------------
553 
554  fait_nphi() ;
555 
556  }
557 
558  //-----------------------------------------
559  // Determination of the fluid velociy U
560  //-----------------------------------------
561 
562  if (mer > mer_fix_omega + delta_mer_kep) {
563 
564  omega *= fact_omega ; // Increase of the angular velocity if
565  } // fact_omega != 1
566 
567  bool omega_trop_grand = false ;
568  bool kepler = true ;
569 
570  while ( kepler ) {
571 
572  // Possible decrease of Omega to ensure a velocity < c
573 
574  bool superlum = true ;
575 
576  while ( superlum ) {
577 
578  // New fluid velocity U :
579 
580  Scalar tmp = omega - nphi ;
581  tmp.annule_domain(nzm1) ;
582  tmp.std_spectral_base() ;
583 
584  tmp.mult_rsint() ; // Multiplication by r sin(theta)
585 
586  uuu = bbb / nn * tmp ;
587 
588  if (uuu.get_etat() == ETATQCQ) {
589  // Same basis as (Omega -N^phi) r sin(theta) :
590  (uuu.set_spectral_va()).set_base( tmp.get_spectral_va().get_base() ) ;
591  }
592 
593  // Is the new velocity larger than c in the equatorial plane ?
594 
595  superlum = false ;
596 
597  for (int l=0; l<nzet; l++) {
598  for (int i=0; i<mg->get_nr(l); i++) {
599 
600  double u1 = uuu.val_grid_point(l, 0, j_b, i) ;
601  if (u1 >= 1.) { // superluminal velocity
602  superlum = true ;
603  cout << "U > c for l, i : " << l << " " << i
604  << " U = " << u1 << endl ;
605  }
606  }
607  }
608  if ( superlum ) {
609  cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ;
610  omega /= fact_omega ; // Decrease of Omega
611  cout << "New rotation frequency : "
612  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
613  omega_trop_grand = true ;
614  }
615  } // end of while ( superlum )
616 
617 
618  // New computation of U (which this time is not superluminal)
619  // as well as of gam_euler, ener_euler, etc...
620  // -----------------------------------
621 
622  hydro_euler() ;
623 
624 
625  //------------------------------------------------------
626  // First integral of motion
627  //------------------------------------------------------
628 
629  // Centrifugal potential :
630  if (relativistic) {
631  mlngamma = - log( gam_euler ) ;
632  }
633  else {
634  mlngamma = - 0.5 * uuu*uuu ;
635  }
636 
637  // Equatorial values of various potentials :
638  double nuf_b = nuf.val_grid_point(l_b, k_b, j_b, i_b) ;
639  double nuq_b = nuq.val_grid_point(l_b, k_b, j_b, i_b) ;
640  double mlngamma_b = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ;
641 
642  // Central values of various potentials :
643  double nuf_c = nuf.val_grid_point(0,0,0,0) ;
644  double nuq_c = nuq.val_grid_point(0,0,0,0) ;
645  double mlngamma_c = 0 ;
646 
647  // Scale factor to ensure that the enthalpy is equal to ent_b at
648  // the equator
649  double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
650  + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
651  alpha_r = sqrt(alpha_r2) ;
652  cout << "alpha_r = " << alpha_r << endl ;
653 
654  // Readjustment of nu :
655  // -------------------
656 
657  logn = alpha_r2 * nuf + nuq ;
658  double nu_c = logn.val_grid_point(0,0,0,0) ;
659 
660  // First integral --> enthalpy in all space
661  //-----------------
662 
663  ent = (ent_c + nu_c + mlngamma_c) - logn - mlngamma ;
664 
665  // Test: is the enthalpy negative somewhere in the equatorial plane
666  // inside the star ? If yes, this means that the Keplerian velocity
667  // has been overstep.
668 
669  kepler = false ;
670  for (int l=0; l<nzet; l++) {
671  int imax = mg->get_nr(l) - 1 ;
672  if (l == l_b) imax-- ; // The surface point is skipped
673  for (int i=0; i<imax; i++) {
674  if ( ent.val_grid_point(l, 0, j_b, i) < 0. ) {
675  kepler = true ;
676  cout << "ent < 0 for l, i : " << l << " " << i
677  << " ent = " << ent.val_grid_point(l, 0, j_b, i) << endl ;
678  }
679  }
680  }
681 
682  if ( kepler ) {
683  cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
684  omega /= fact_omega ; // Omega is decreased
685  cout << "New rotation frequency : "
686  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
687  omega_trop_grand = true ;
688  if (fact_omega == 1.) throw true ;
689  }
690 
691  } // End of while ( kepler )
692 
693  if ( omega_trop_grand ) { // fact_omega is decreased for the
694  // next step
695  fact_omega = sqrt( fact_omega ) ;
696  cout << "**** New fact_omega : " << fact_omega << endl ;
697  }
698 
699  //----------------------------------------------------
700  // Adaptation of the mapping to the new enthalpy field
701  //----------------------------------------------------
702 
703  // Shall the adaptation be performed (cusp) ?
704  // ------------------------------------------
705 
706  double dent_eq = ent.dsdr().val_grid_point(l_b, k_b, j_b, i_b) ;
707  double dent_pole = ent.dsdr().val_grid_point(l_b, k_b, 0, i_b) ;
708  double rap_dent = fabs( dent_eq / dent_pole ) ;
709  cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
710 
711  if ( rap_dent < thres_adapt ) {
712  adapt_flag = 0 ; // No adaptation of the mapping
713  cout << "******* FROZEN MAPPING *********" << endl ;
714  }
715  else{
716  adapt_flag = 1 ; // The adaptation of the mapping is to be
717  // performed
718  }
719 
720  mp_prev = mp_et ;
721 
722  Cmp cent(ent) ;
723  mp.adapt(cent, par_adapt) ;
724 
725  //----------------------------------------------------
726  // Computation of the enthalpy at the new grid points
727  //----------------------------------------------------
728 
729  mp_prev.homothetie(alpha_r) ;
730 
731  mp.reevaluate(&mp_prev, nzet+1, cent) ;
732  ent = cent ;
733 
734  //----------------------------------------------------
735  // Equation of state
736  //----------------------------------------------------
737 
738  equation_of_state() ; // computes new values for nbar (n), ener (e)
739  // and press (p) from the new ent (H)
740 
741  //---------------------------------------------------------
742  // Matter source terms in the gravitational field equations
743  //---------------------------------------------------------
744 
745  //## Computation of tnphi and nphi from the Cartesian components
746  // of the shift for the test in hydro_euler():
747 
748  fait_nphi() ;
749 
750  hydro_euler() ; // computes new values for ener_euler (E),
751  // s_euler (S) and u_euler (U^i)
752 
753  if (relativistic) {
754 
755  //-------------------------------------------------------
756  // 2-D Poisson equation for tggg
757  //-------------------------------------------------------
758 
759  Cmp csource_tggg(source_tggg) ;
760  Cmp ctggg(tggg) ;
761  mp.poisson2d(csource_tggg, mp.cmp_zero(), par_poisson_tggg,
762  ctggg) ;
763  tggg = ctggg ;
764 
765 
766  //-------------------------------------------------------
767  // 2-D Poisson equation for dzeta
768  //-------------------------------------------------------
769 
770  Cmp csource_dzf(source_dzf) ;
771  Cmp csource_dzq(source_dzq) ;
772  Cmp cdzeta(dzeta) ;
773  mp.poisson2d(csource_dzf, csource_dzq, par_poisson_dzeta,
774  cdzeta) ;
775  dzeta = cdzeta ;
776 
777  err_grv2 = lbda_grv2 - 1;
778  cout << "GRV2: " << err_grv2 << endl ;
779 
780  }
781  else {
782  err_grv2 = grv2() ;
783  }
784 
785 
786  //---------------------------------------
787  // Computation of the metric coefficients (except for N^phi)
788  //---------------------------------------
789 
790  // Relaxations on nu and dzeta :
791 
792  if (mer >= 10) {
793  logn = relax * logn + relax_prev * logn_prev ;
794 
795  dzeta = relax * dzeta + relax_prev * dzeta_prev ;
796  }
797 
798  // Update of the metric coefficients N, A, B and computation of K_ij :
799 
800  update_metric() ;
801 
802  //-----------------------
803  // Informations display
804  //-----------------------
805 
806  partial_display(cout) ;
807  fichfreq << " " << setprecision(16) << omega / (2*M_PI) * f_unit ;
808  fichevol << " " << setprecision(16) << rap_dent ;
809  fichevol << " " << setprecision(16) << ray_pole() / ray_eq() ;
810  fichevol << " " << setprecision(16) << ent_c ;
811 
812  //-----------------------------------------
813  // Convergence towards a given baryon mass
814  //-----------------------------------------
815 
816  if (mer > mer_mass) {
817 
818  double xx ;
819  if (mbar_wanted > 0.) {
820  xx = mass_b() / mbar_wanted - 1. ;
821  cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx
822  << endl ;
823  }
824  else{
825  xx = mass_g() / fabs(mbar_wanted) - 1. ;
826  cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx
827  << endl ;
828  }
829  double xprog = ( mer > 2*mer_mass) ? 1. :
830  double(mer-mer_mass)/double(mer_mass) ;
831  xx *= xprog ;
832  double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
833  double fact = pow(ax, aexp_mass) ;
834  cout << " xprog, xx, ax, fact : " << xprog << " " <<
835  xx << " " << ax << " " << fact << endl ;
836 
837  if ( change_ent ) {
838  ent_c *= fact ;
839  }
840  else {
841  if (mer%4 == 0) omega *= fact ;
842  }
843  }
844 
845 
846  //------------------------------------------------------------
847  // Relative change in enthalpy with respect to previous step
848  //------------------------------------------------------------
849 
850  Tbl diff_ent_tbl = diffrel( ent, ent_prev ) ;
851  diff_ent = diff_ent_tbl(0) ;
852  for (int l=1; l<nzet; l++) {
853  diff_ent += diff_ent_tbl(l) ;
854  }
855  diff_ent /= nzet ;
856 
857  fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
858  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
859  fichconv << " " << log10( fabs(max_triax) + 1.e-16 ) ;
860 
861  vit_triax = 0 ;
862  if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
863  vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
864  }
865 
866  fichconv << " " << vit_triax ;
867 
868  //------------------------------
869  // Recycling for the next step
870  //------------------------------
871 
872  ent_prev = ent ;
873  logn_prev = logn ;
874  dzeta_prev = dzeta ;
875  max_triax_prev = max_triax ;
876 
877  fichconv << endl ;
878  fichfreq << endl ;
879  fichevol << endl ;
880  fichconv.flush() ;
881  fichfreq.flush() ;
882  fichevol.flush() ;
883 
884  } // End of main loop
885 
886  //=========================================================================
887  // End of iteration
888  //=========================================================================
889 
890  ssjm1_tggg = cssjm1_tggg ;
891  ssjm1_khi = cssjm1_khi ;
892  for (int i=1; i<=3; i++) {
893  ssjm1_wshift.set(i) = cssjm1_wshift(i-1) ;
894  }
895 
896  fichconv.close() ;
897  fichfreq.close() ;
898  fichevol.close() ;
899 
900 }
901 }
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.
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
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
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 hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
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
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
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Definition: star_rot.h:220
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
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
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
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
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
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