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