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