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