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