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