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