LORENE
et_bfrot_equilibre.C
1/*
2 * Method of class Et_rot_bifluid to compute a static spherical configuration.
3 *
4 * (see file etoile.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2001 Jerome Novak
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
30char et_bfrot_equilibre_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_equilibre.C,v 1.20 2015/06/10 14:39:17 a_sourie Exp $" ;
31
32/*
33 * $Id: et_bfrot_equilibre.C,v 1.20 2015/06/10 14:39:17 a_sourie Exp $
34 * $Log: et_bfrot_equilibre.C,v $
35 * Revision 1.20 2015/06/10 14:39:17 a_sourie
36 * New class Eos_bf_tabul for tabulated 2-fluid EoSs and associated functions for the computation of rotating stars with such EoSs.
37 *
38 * Revision 1.19 2014/10/13 08:52:54 j_novak
39 * Lorene classes and functions now belong to the namespace Lorene.
40 *
41 * Revision 1.18 2014/10/06 15:13:07 j_novak
42 * Modified #include directives to use c++ syntax.
43 *
44 * Revision 1.17 2006/03/13 10:02:27 j_novak
45 * Added things for triaxial perturbations.
46 *
47 * Revision 1.16 2004/09/01 10:56:05 r_prix
48 * added option of converging baryon-mass to equilibrium_bi()
49 *
50 * Revision 1.15 2004/08/30 09:54:20 r_prix
51 * experimental version of Kepler-limit finder for 2-fluid stars
52 *
53 * Revision 1.14 2004/03/25 10:29:03 j_novak
54 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
55 *
56 * Revision 1.13 2003/12/11 12:43:35 r_prix
57 * activated adaptive grid for 2-fluid star (taken from Etoile_rot)
58 *
59 * Revision 1.12 2003/12/04 14:28:26 r_prix
60 * allow for the case of "slow-rot-style" EOS inversion, in which we need to adapt
61 * the inner domain to n_outer=0 instead of mu_outer=0 ...
62 * (this should only be used for comparison to analytic slow-rot solution!)
63 *
64 * Revision 1.11 2003/11/25 12:49:44 j_novak
65 * Modified headers to compile on IRIX. Changed Mapping to be Map_af (speed
66 * enhancement).
67 *
68 * Revision 1.10 2003/11/20 14:01:25 r_prix
69 * changed member names to better conform to Lorene coding standards:
70 * J_euler -> j_euler, EpS_euler -> enerps_euler, Delta_car -> delta_car
71 *
72 * Revision 1.9 2003/11/19 22:01:57 e_gourgoulhon
73 * -- Relaxation on logn and dzeta performed only if mer >= 10.
74 * -- err_grv2 is now evaluated also in the Newtonian case.
75 *
76 * Revision 1.8 2003/11/18 18:38:11 r_prix
77 * use of new member EpS_euler: matter sources in equilibrium() and global quantities
78 * no longer distinguish Newtonian/relativistic, as all terms should have the right limit...
79 *
80 * Revision 1.7 2003/11/17 13:49:43 r_prix
81 * - moved superluminal check into hydro_euler()
82 * - removed some warnings
83 *
84 * Revision 1.6 2003/11/13 12:14:35 r_prix
85 * *) removed all use of etoile-type specific u_euler, press
86 * and use 3+1 components of Tmunu instead
87 *
88 * Revision 1.5 2002/10/16 14:36:35 j_novak
89 * Reorganization of #include instructions of standard C++, in order to
90 * use experimental version 3 of gcc.
91 *
92 * Revision 1.4 2002/04/05 09:09:36 j_novak
93 * The inversion of the EOS for 2-fluids polytrope has been modified.
94 * Some errors in the determination of the surface were corrected.
95 *
96 * Revision 1.3 2002/01/16 15:03:28 j_novak
97 * *** empty log message ***
98 *
99 * Revision 1.2 2002/01/03 15:30:28 j_novak
100 * Some comments modified.
101 *
102 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
103 * LORENE
104 *
105 * Revision 1.1 2001/06/22 15:40:06 novak
106 * Initial revision
107 *
108 *
109 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_equilibre.C,v 1.20 2015/06/10 14:39:17 a_sourie Exp $
110 *
111 */
112
113// Headers C
114#include <cmath>
115
116// Headers Lorene
117#include "et_rot_bifluid.h"
118#include "param.h"
119#include "unites.h"
120
121#include "graphique.h"
122#include "utilitaires.h"
123
124namespace Lorene {
125
126//-------------------------------------------------------------------------
127//-------------------------------------------------------------------------
128
129// Axial Equilibrium
130
131//-------------------------------------------------------------------------
132//-------------------------------------------------------------------------
133
135(double ent_c, double ent2_c, double omega0, double omega20,
136 const Tbl& ent_limit, const Tbl& ent2_limit, const Itbl& icontrol,
137 const Tbl& control, Tbl& diff,
138 int mer_mass, double mbar1_wanted, double mbar2_wanted, double aexp_mass)
139{
140
141 // Fundamental constants and units
142 // -------------------------------
143 using namespace Unites ;
144
145 // For the display
146 // ---------------
147 char display_bold[]="x[1m" ; display_bold[0] = 27 ;
148 char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
149
150 // Grid parameters
151 // ---------------
152
153 const Mg3d* mg = mp.get_mg() ;
154 int nz = mg->get_nzone() ; // total number of domains
155
156 // The following is required to initialize mp_prev as a Map_af
157 Map_et& mp_et = dynamic_cast<Map_et&>(mp) ; // reference
158
159 // Index of the point at phi=0, theta=pi/2 at the surface of the star:
160 assert(mg->get_type_t() == SYM) ;
161 int l_b = nzet - 1 ;
162 int i_b = mg->get_nr(l_b) - 1 ;
163 int j_b = mg->get_nt(l_b) - 1 ;
164 int k_b = 0 ;
165
166 // Value of the enthalpies defining the surface of each fluid
167 double ent1_b = ent_limit(nzet-1) ;
168 double ent2_b = ent2_limit(nzet-1) ;
169
170 // This value is chosen so that the grid contain both fluids
171// double ent_b = (ent1_b > ent2_b ? ent1_b : ent2_b) ;
172
173 // Parameters to control the iteration
174 // -----------------------------------
175
176 int mer_max = icontrol(0) ;
177 int mer_rot = icontrol(1) ;
178 int mer_change_omega = icontrol(2) ;
179 int mer_fix_omega = icontrol(3) ;
180 int mermax_poisson = icontrol(4) ;
181 int nzadapt = icontrol(5); // number of domains for adaptive grid
182 int kepler_fluid = icontrol(6); // fluid-index for kepler-limit (0=none,3=both)
183 int kepler_wait_steps = icontrol(7);
184 int mer_triax = icontrol(8) ;
185
186 int niter ;
187
188 // Protections:
189 if (mer_change_omega < mer_rot) {
190 cout << "Et_rot_bifluid::equilibrium: mer_change_omega < mer_rot !" << endl ;
191 cout << " mer_change_omega = " << mer_change_omega << endl ;
192 cout << " mer_rot = " << mer_rot << endl ;
193 abort() ;
194 }
195 if (mer_fix_omega < mer_change_omega) {
196 cout << "Et_rot_bifluid::equilibrium: mer_fix_omega < mer_change_omega !"
197 << endl ;
198 cout << " mer_fix_omega = " << mer_fix_omega << endl ;
199 cout << " mer_change_omega = " << mer_change_omega << endl ;
200 abort() ;
201 }
202
203 double precis = control(0) ;
204 double omega_ini = control(1) ;
205 double omega2_ini = control(2) ;
206 double relax = control(3) ;
207 double relax_prev = double(1) - relax ;
208 double relax_poisson = control(4) ;
209 // some additional stuff for adaptive grid:
210 double thres_adapt = control(5) ;
211 double precis_adapt = control(6) ;
212 double kepler_factor = control(7);
213 if (kepler_factor <= 1.0)
214 {
215 cout << "ERROR: Kepler factor has to be greater than 1!!\n";
216 abort();
217 }
218 double ampli_triax = control(8) ;
219
220 // Error indicators
221 // ----------------
222
223 diff.set_etat_qcq() ;
224 double diff_ent ;
225 double& diff_ent1 = diff.set(0) ;
226 double& diff_ent2 = diff.set(1) ;
227 double& diff_nuf = diff.set(2) ;
228 double& diff_nuq = diff.set(3) ;
229 // double& diff_dzeta = diff.set(4) ;
230 // double& diff_ggg = diff.set(5) ;
231 double& diff_shift_x = diff.set(6) ;
232 double& diff_shift_y = diff.set(7) ;
233 double& vit_triax = diff.set(8) ;
234
235 // Parameters for the function Map_et::adapt
236 // -----------------------------------------
237
238 Param par_adapt ;
239 int nitermax = 100 ;
240 int adapt_flag = 1 ; // 1 = performs the full computation,
241 // 0 = performs only the rescaling by
242 // the factor alpha_r
243 int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
244 // isosurfaces
245 double alpha_r ;
246 double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
247
248 par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
249 // locate zeros by the secant method
250 par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment
251 // to the isosurfaces of ent is to be
252 // performed
253 par_adapt.add_int(nz_search, 2) ; // number of domains to search for
254 // the enthalpy isosurface
255 par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
256 // 0 = performs only the rescaling by
257 // the factor alpha_r
258 par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
259 // (theta_*, phi_*)
260 par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
261 // (theta_*, phi_*)
262
263 par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
264 // the secant method
265
266 par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in
267 // the determination of zeros by
268 // the secant method
269 par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping, 0 = contracting mapping
270
271 par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
272 // distances will be multiplied
273
274 par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
275 // to define the isosurfaces.
276
277
278 // Parameters for the function Map_et::poisson for nuf
279 // ----------------------------------------------------
280
281 double precis_poisson = 1.e-16 ;
282
283 Param par_poisson_nuf ;
284 par_poisson_nuf.add_int(mermax_poisson, 0) ; // maximum number of iterations
285 par_poisson_nuf.add_double(relax_poisson, 0) ; // relaxation parameter
286 par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
287 par_poisson_nuf.add_int_mod(niter, 0) ; // number of iterations actually used
288 par_poisson_nuf.add_cmp_mod( ssjm1_nuf ) ;
289
290 Param par_poisson_nuq ;
291 par_poisson_nuq.add_int(mermax_poisson, 0) ; // maximum number of iterations
292 par_poisson_nuq.add_double(relax_poisson, 0) ; // relaxation parameter
293 par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
294 par_poisson_nuq.add_int_mod(niter, 0) ; // number of iterations actually used
295 par_poisson_nuq.add_cmp_mod( ssjm1_nuq ) ;
296
297 Param par_poisson_tggg ;
298 par_poisson_tggg.add_int(mermax_poisson, 0) ; // maximum number of iterations
299 par_poisson_tggg.add_double(relax_poisson, 0) ; // relaxation parameter
300 par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
301 par_poisson_tggg.add_int_mod(niter, 0) ; // number of iterations actually used
302 par_poisson_tggg.add_cmp_mod( ssjm1_tggg ) ;
303 double lambda_tggg ;
304 par_poisson_tggg.add_double_mod( lambda_tggg ) ;
305
306 Param par_poisson_dzeta ;
307 double lbda_grv2 ;
308 par_poisson_dzeta.add_double_mod( lbda_grv2 ) ;
309
310 // Parameters for the function Tenseur::poisson_vect
311 // -------------------------------------------------
312
313 Param par_poisson_vect ;
314
315 par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
316 par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
317 par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
318 par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
319 par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
320 par_poisson_vect.add_int_mod(niter, 0) ;
321
322
323 // Initializations
324 // ---------------
325
326 // Initial angular velocities
327 omega = 0 ;
328 omega2 = 0 ;
329
330 double accrois_omega = (omega0 - omega_ini) /
331 double(mer_fix_omega - mer_change_omega) ;
332 double accrois_omega2 = (omega20 - omega2_ini) /
333 double(mer_fix_omega - mer_change_omega) ;
334
335
336 update_metric() ; // update of the metric coefficients
337
338 equation_of_state() ; // update of the densities, pressure, etc...
339
340 hydro_euler() ; // update of the hydro quantities relative to the
341 // Eulerian observer
342
343 // Quantities at the previous step :
344
345 Map_et mp_prev = mp_et;
346 Tenseur ent_prev = ent ;
347 Tenseur ent2_prev = ent2 ;
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 Tenseur mlngamma2(mp) ; // centrifugal potential
361
362 Tenseur *outer_ent_p; // pointer to the enthalpy field of the outer fluid
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) {
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 << "# f1 [Hz] f2 [Hz]" << endl ;
393
394 ofstream fichevol("evolution.d") ; // Output file for various quantities
395 fichevol << "# r_pole/r_eq ent_c ent2_c" << endl ;
396
397 diff_ent = 1 ;
398 double err_grv2 = 1 ;
399 double max_triax_prev = 0 ; // Triaxial amplitude at previous step
400
401 //=========================================================================
402 // Start of iteration
403 //=========================================================================
404
405 for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
406
407 cout << "-----------------------------------------------" << endl ;
408 cout << "step: " << mer << endl ;
409 cout << "diff_ent = " << display_bold << diff_ent << display_normal
410 << endl ;
411 cout << "err_grv2 = " << err_grv2 << endl ;
412 fichconv << mer ;
413 fichfreq << mer ;
414 fichevol << mer ;
415
416 if (mer >= mer_rot) {
417
418 if (mer < mer_change_omega) {
419 omega = omega_ini ;
420 omega2 = omega2_ini ;
421 }
422 else {
423 if (mer <= mer_fix_omega) {
424 omega = omega_ini + accrois_omega *
425 (mer - mer_change_omega) ;
426 omega2 = omega2_ini + accrois_omega2 *
427 (mer - mer_change_omega) ;
428 }
429 }
430
431 }
432
433 //-----------------------------------------------
434 // Sources of the Poisson equations
435 //-----------------------------------------------
436
437 // Source for nu
438 // -------------
439 Tenseur beta = log(bbb) ;
440 beta.set_std_base() ;
441
442 // common source term for relativistic and Newtonian ! (enerps_euler has the right limit)
443 source_nuf = qpig * a_car * enerps_euler;
444
445 if (relativistic)
447 else
448 source_nuq = 0 ;
449
450 source_nuf.set_std_base() ;
451 source_nuq.set_std_base() ;
452
453 if (relativistic) {
454 // Source for dzeta
455 // ----------------
456 source_dzf = 2 * qpig * a_car * sphph_euler;
457 source_dzf.set_std_base() ;
458
459 source_dzq = 1.5 * ak_car - flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() ) ;
460 source_dzq.set_std_base() ;
461
462 // Source for tggg
463 // ---------------
464
465 source_tggg = 2 * qpig * nnn * a_car * bbb * (s_euler - sphph_euler);
466 source_tggg.set_std_base() ;
467
468 (source_tggg.set()).mult_rsint() ;
469
470
471 // Source for shift
472 // ----------------
473
474 // Matter term
475 source_shift = (-4*qpig) * nnn * a_car * j_euler;
476
477 // Quadratic terms:
478 Tenseur vtmp = 3 * beta.gradient_spher() - logn.gradient_spher() ;
480
481 Tenseur squad = 2 * 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 diff_nuf = 0 ;
505
506 //---------------------------------------
507 // Triaxial perturbation of nuf
508 //---------------------------------------
509
510 if (mer == mer_triax) {
511
512 if ( mg->get_np(0) == 1 ) {
513 cout <<
514 "Et_rot_bifluid::equilibrium: np must be stricly greater than 1"
515 << endl << " to set a triaxial perturbation !" << endl ;
516 abort() ;
517 }
518
519 const Coord& phi = mp.phi ;
520 const Coord& sint = mp.sint ;
521 Cmp perturb(mp) ;
522 perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
523 nuf.set() = nuf() * perturb ;
524
525 nuf.set_std_base() ; // set the bases for spectral expansions
526 // to be the standard ones for a
527 // scalar field
528
529 }
530
531 // Monitoring of the triaxial perturbation
532 // ---------------------------------------
533
534 Valeur& va_nuf = nuf.set().va ;
535 va_nuf.coef() ; // Computes the spectral coefficients
536 double max_triax = 0 ;
537
538 if ( mg->get_np(0) > 1 ) {
539
540 for (int l=0; l<nz; l++) { // loop on the domains
541 for (int j=0; j<mg->get_nt(l); j++) {
542 for (int i=0; i<mg->get_nr(l); i++) {
543
544 // Coefficient of cos(2 phi) :
545 double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ;
546
547 // Coefficient of sin(2 phi) :
548 double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ;
549
550 double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
551
552 max_triax = ( xx > max_triax ) ? xx : max_triax ;
553 }
554 }
555 }
556 }
557 cout << "Triaxial part of nuf : " << max_triax << endl ;
558
559 if (relativistic) {
560
561 //----------------------------------------------
562 // Resolution of the Poisson equation for nuq
563 //----------------------------------------------
564
565 source_nuq().poisson(par_poisson_nuq, nuq.set()) ;
566
567// cout << "Test of the Poisson equation for nuq :" << endl ;
568// err = source_nuq().test_poisson(nuq(), cout, true) ;
569// diff_nuq = err(0, 0) ;
570 diff_nuq = 0 ;
571
572 //---------------------------------------------------------
573 // Resolution of the vector Poisson equation for the shift
574 //---------------------------------------------------------
575
576 if (source_shift.get_etat() != ETATZERO) {
577
578 for (int i=0; i<3; i++) {
579 if(source_shift(i).dz_nonzero()) {
580 assert( source_shift(i).get_dzpuis() == 4 ) ;
581 }
582 else{
583 (source_shift.set(i)).set_dzpuis(4) ;
584 }
585 }
586
587 }
588
589 double lambda_shift = double(1) / double(3) ;
590
591 if ( mg->get_np(0) == 1 ) {
592 lambda_shift = 0 ;
593 }
594
595 source_shift.poisson_vect(lambda_shift, par_poisson_vect,
597
598// cout << "Test of the Poisson equation for shift_x :" << endl ;
599// err = source_shift(0).test_poisson(shift(0), cout, true) ;
600// diff_shift_x = err(0, 0) ;
601 diff_shift_x = 0 ;
602
603// cout << "Test of the Poisson equation for shift_y :" << endl ;
604// err = source_shift(1).test_poisson(shift(1), cout, true) ;
605// diff_shift_y = err(0, 0) ;
606 diff_shift_y = 0 ;
607
608 // Computation of tnphi and nphi from the Cartesian components
609 // of the shift
610 // -----------------------------------------------------------
611
612 fait_nphi() ;
613
614 }
615
616
617 //----------------------------------------
618 // Shall we search for the Kepler limit?
619 //----------------------------------------
620 bool kepler = false;
621 bool too_fast = false;
622
623 if ( (kepler_fluid > 0) && (mer > mer_fix_omega + kepler_wait_steps) )
624 {
625 if (kepler_fluid & 0x01)
626 omega *= kepler_factor;
627 if (kepler_fluid & 0x02)
628 omega2 *= kepler_factor;
629 }
630
631
632 // ============================================================
633 kepler = true;
634 while (kepler)
635 {
636
637 // New computation of delta_car, gam_euler, enerps_euler etc...
638 // ------------------------------------------------------
639 hydro_euler() ;
640
641
642 //------------------------------------------------------
643 // First integral of motion
644 //------------------------------------------------------
645
646 // Centrifugal potential :
647 if (relativistic) {
648 mlngamma = - log( gam_euler ) ;
649 mlngamma2 = - log( gam_euler2) ;
650 }
651 else {
652 mlngamma = - 0.5 * uuu*uuu ;
653 mlngamma2 = -0.5 * uuu2*uuu2 ;
654 }
655
656 // Central values of various potentials :
657 double nuf_c = nuf()(0,0,0,0) ;
658 double nuq_c = nuq()(0,0,0,0) ;
659
660 // Scale factor to ensure that the enthalpy is equal to ent_b at
661 // the equator for the "outer" fluid
662 double alpha_r2 = 0;
663
664 int j=j_b;
665
666 // Boundary values of various potentials :
667 double nuf_b = nuf()(l_b, k_b, j, i_b) ;
668 double nuq_b = nuq()(l_b, k_b, j, i_b) ;
669 double mlngamma_b = mlngamma()(l_b, k_b, j, i_b) ;
670 double mlngamma2_b = mlngamma2()(l_b, k_b, j, i_b) ;
671
672
673 // RP: "hack": adapt the radius correctly if using "slow-rot-style" EOS inversion
674 //
675 if ( eos.identify() == 2 ) // only applies to Eos_bf_poly_newt
676 {
677 const Eos_bf_poly_newt &eos0 = dynamic_cast<const Eos_bf_poly_newt&>(eos);
678 if (eos0.get_typeos() == 5)
679 {
680 double vn_b = uuu()(l_b, k_b, j, i_b);
681 double vp_b = uuu2()(l_b, k_b, j, i_b);
682 double D2_b = (vp_b - vn_b)*(vp_b - vn_b);
683 double kdet = eos0.get_kap3() + eos0.get_beta()*D2_b;
684 double kaps1 = kdet / ( eos0.get_kap2() - kdet );
685 double kaps2 = kdet / ( eos0.get_kap1() - kdet );
686
687 ent1_b = kaps1 * ( ent2_c - ent_c - mlngamma2_b + mlngamma_b );
688 ent2_b = kaps2 * ( ent_c - ent2_c - mlngamma_b + mlngamma2_b );
689
690 cout << "**********************************************************************\n";
691 cout << "DEBUG: Rescaling domain for slow-rot-style EOS inversion \n";
692 cout << "DEBUG: ent1_b = " << ent1_b << "; ent2_b = " << ent2_b << endl;
693 cout << "**********************************************************************\n";
694
695 adapt_flag = 0; // don't do adaptive-grid if using slow-rot-style inversion!
696 }
697 }
698
699 double alpha1_r2 = ( ent_c - ent1_b - mlngamma_b + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
700 double alpha2_r2 = ( ent2_c - ent2_b - mlngamma2_b + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
701
702 cout << "DEBUG: j= "<< j<<" ; alpha1 = " << alpha1_r2 <<" ; alpha2 = " << alpha2_r2 << endl;
703
704 int outer_fluid = (alpha1_r2 > alpha2_r2) ? 1 : 2; // index of 'outer' fluid (at equator!)
705
706 outer_ent_p = (outer_fluid == 1) ? (&ent) : (&ent2);
707
708 alpha_r2 = (outer_fluid == 1) ? alpha1_r2 : alpha2_r2 ;
709
710 alpha_r = sqrt(alpha_r2);
711
712 cout << "alpha_r = " << alpha_r << endl ;
713
714 // Readjustment of nu :
715 // -------------------
716
717 logn = alpha_r2 * nuf + nuq ;
718 double nu_c = logn()(0,0,0,0) ;
719
720
721 // First integral --> enthalpy in all space
722 //-----------------
723
724 ent = (ent_c + nu_c) - logn - mlngamma ;
725 ent2 = (ent2_c + nu_c) - logn - mlngamma2 ;
726
727
728 // now let's try to figure out if we have overstepped the Kepler-limit
729 // (FIXME) we assume that the enthalpy of the _outer_ fluid being negative
730 // inside the star
731 kepler = false;
732 for (int l=0; l<nzet; l++) {
733 int imax = mg->get_nr(l) - 1 ;
734 if (l == l_b) imax-- ; // The surface point is skipped
735 for (int i=0; i<imax; i++) {
736 if ( (*outer_ent_p)()(l, 0, j_b, i) < 0. ) {
737 kepler = true;
738 cout << "(outer) ent < 0 for l, i : " << l << " " << i
739 << " ent = " << (*outer_ent_p)()(l, 0, j_b, i) << endl ;
740 }
741 }
742 }
743
744 if ( kepler )
745 {
746 cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
747 if (kepler_fluid & 0x01)
748 omega /= kepler_factor ; // Omega is decreased
749 if (kepler_fluid & 0x02)
750 omega2 /= kepler_factor;
751
752 cout << "New rotation frequencies : "
753 << "Omega = " << omega/(2.*M_PI) * f_unit << " Hz; "
754 << "Omega2 = " << omega2/(2.*M_PI) * f_unit << endl ;
755
756 too_fast = true;
757 }
758
759 } /* while kepler */
760
761
762 if ( too_fast )
763 { // fact_omega is decreased for the next step
764 kepler_factor = sqrt( kepler_factor ) ;
765 cout << "**** New fact_omega : " << kepler_factor << endl ;
766 }
767 // ============================================================
768
769
770 // Cusp-check: shall the adaptation still be performed?
771 // ------------------------------------------
772 double dent_eq = (*outer_ent_p)().dsdr()(l_b, k_b, j_b, i_b) ;
773 double dent_pole = (*outer_ent_p)().dsdr()(l_b, k_b, 0, i_b) ;
774 double rap_dent = fabs( dent_eq / dent_pole ) ;
775 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
776
777 if ( rap_dent < thres_adapt ) {
778 adapt_flag = 0 ; // No adaptation of the mapping
779 cout << "******* FROZEN MAPPING *********" << endl ;
780 }
781
782 // Rescaling of the grid and adaption to (outer) enthalpy surface
783 //---------------------------------------
784 if (adapt_flag && (nzadapt > 0) )
785 {
786 mp_prev = mp_et ;
787
788 mp.adapt( (*outer_ent_p)(), par_adapt) ;
789
790 mp_prev.homothetie(alpha_r) ;
791
792 mp.reevaluate(&mp_prev, nzet+1, ent.set()) ;
793 mp.reevaluate(&mp_prev, nzet+1, ent2.set()) ;
794 }
795 else
796 mp.homothetie (alpha_r);
797
798
799 //----------------------------------------------------
800 // Equation of state
801 //----------------------------------------------------
802
803 equation_of_state() ; // computes new values for nbar1,2 , ener (e)
804 // and press (p) from the new ent,ent2
805
806 //---------------------------------------------------------
807 // Matter source terms in the gravitational field equations
808 //---------------------------------------------------------
809
810 //## Computation of tnphi and nphi from the Cartesian components
811 // of the shift for the test in hydro_euler():
812
813 fait_nphi() ;
814
815 hydro_euler() ; // computes new values for ener_euler (E),
816 // s_euler (S) and u_euler (U^i)
817
818 if (relativistic) {
819
820 //-------------------------------------------------------
821 // 2-D Poisson equation for tggg
822 //-------------------------------------------------------
823
824 mp.poisson2d(source_tggg(), mp.cmp_zero(), par_poisson_tggg, tggg.set()) ;
825
826 //-------------------------------------------------------
827 // 2-D Poisson equation for dzeta
828 //-------------------------------------------------------
829
830 mp.poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta, dzeta.set()) ;
831
832 err_grv2 = lbda_grv2 - 1;
833 cout << "GRV2: " << err_grv2 << endl ;
834
835 }
836 else {
837 err_grv2 = grv2() ;
838 }
839
840
841 //---------------------------------------
842 // Computation of the metric coefficients (except for N^phi)
843 //---------------------------------------
844
845 // Relaxations on nu and dzeta :
846
847 if (mer >= 10) {
848 logn = relax * logn + relax_prev * logn_prev ;
849
850 dzeta = relax * dzeta + relax_prev * dzeta_prev ;
851 }
852
853 // Update of the metric coefficients N, A, B and computation of K_ij :
854
855 update_metric() ;
856
857 //-----------------------
858 // Informations display
859 //-----------------------
860
861 // partial_display(cout) ;
862 fichfreq << " " << omega / (2*M_PI) * f_unit ;
863 fichfreq << " " << omega2 / (2*M_PI) * f_unit ;
864 fichevol << " " << ray_pole() / ray_eq() ;
865 fichevol << " " << ent_c ;
866 fichevol << " " << ent2_c ;
867
868
869 //-----------------------------------------
870 // Convergence towards given baryon masses (if mer_mass > 0)
871 //-----------------------------------------
872
873 // If we want to impose baryonic masses for both fluids.
874 //Be careful, the code acts on mu_n and mu_p (at the center)
875 // -> beta equilibrium can be not verified
876 cout << "DEBUG MODE : mbar1_wanted : " << mbar1_wanted << endl ;
877 cout << "DEBUG MODE : mbar2_wanted : " << mbar2_wanted << endl ;
878 if (mbar2_wanted != 0 )
879
880 {
881 if ((mer_mass>0) && (mer > mer_mass)) {
882
883 double xx, xprog, ax, fact;
884
885 // fluid 1
886 xx = mass_b1() / mbar1_wanted - 1. ;
887 cout << "Discrep. baryon mass1 <-> wanted bar. mass1 : " << xx << endl ;
888
889 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
890 xx *= xprog ;
891 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
892 fact = pow(ax, aexp_mass) ;
893 cout << "Fluid1: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
894 ent_c *= fact ;
895
896 // fluid 2
897 xx = mass_b2() / mbar2_wanted - 1. ;
898 cout << "Discrep. baryon mass2 <-> wanted bar. mass2 : " << xx << endl ;
899
900 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
901 xx *= xprog ;
902 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
903 fact = pow(ax, aexp_mass) ;
904 cout << "Fluid2: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
905 ent2_c *= fact ;
906 cout << "H1c = " << ent_c << " H2c = " << ent2_c << endl ;
907 }
908 }
909
910 else {
911 // If we want to impose Mb_tot and beta equilibrium
912 // In this case : mbar1_wanted = total baryonic mass wanted
913 // mbar2_wanted must be set to -1 and is not used.
914
915 if ((mer_mass>0) && (mer > mer_mass)) {
916
917 double xx, xprog, ax, fact;
918
919 // total mass
920 xx = mass_b() / mbar1_wanted - 1. ; // mbar1_wanted = " mbar_wanted"
921 cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx << endl ;
922
923 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
924 xx *= xprog ;
925 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
926 fact = pow(ax, aexp_mass) ;
927 cout << "Fluid1: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
928 ent_c *= fact ;
929
930 double m1 = 1.009000285 ;
931 double m2 = 1.008160139 ;
932 ent2_c = ent_c + log(m1/m2); // to ensure beta_equilibrium
933 cout << "DEBUG MODE : ent_c " << ent_c << endl ;
934 cout << "DEBUG MODE : ent2_c " << ent2_c << endl ;
935 cout << "H1c = " << ent_c << " H2c = " << ent2_c << endl ;
936
937 }
938
939 } /* if mer > mer_mass */
940
941
942
943
944
945 //-------------------------------------------------------------
946 // Relative change in enthalpies with respect to previous step
947 //-------------------------------------------------------------
948
949 Tbl diff_ent_tbl = diffrel( ent(), ent_prev() ) ;
950 diff_ent1 = diff_ent_tbl(0) ;
951 for (int l=1; l<nzet; l++) {
952 diff_ent1 += diff_ent_tbl(l) ;
953 }
954 diff_ent1 /= nzet ;
955 diff_ent_tbl = diffrel( ent2(), ent2_prev() ) ;
956 diff_ent2 = diff_ent_tbl(0) ;
957 for (int l=1; l<nzet; l++) {
958 diff_ent2 += diff_ent_tbl(l) ;
959 }
960 diff_ent2 /= nzet ;
961 diff_ent = 0.5*(diff_ent1 + diff_ent2) ;
962
963 fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
964 fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
965 fichconv << " " << log10( fabs(max_triax) + 1.e-16 ) ;
966
967 vit_triax = 0 ;
968 if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
969 vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
970 }
971
972 fichconv << " " << vit_triax ;
973
974 //------------------------------
975 // Recycling for the next step
976 //------------------------------
977
978 ent_prev = ent ;
979 ent2_prev = ent2 ;
980 logn_prev = logn ;
981 dzeta_prev = dzeta ;
982 max_triax_prev = max_triax ;
983
984 fichconv << endl ;
985 fichfreq << endl ;
986 fichevol << endl ;
987 fichconv.flush() ;
988 fichfreq.flush() ;
989 fichevol.flush() ;
990
991 } // End of main loop
992
993 //=========================================================================
994 // End of iteration
995 //=========================================================================
996
997 fichconv.close() ;
998 fichfreq.close() ;
999 fichevol.close() ;
1000
1001
1002}
1003
1004}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Valeur va
The numerical value of the Cmp
Definition cmp.h:464
Active physical coordinates and mapping derivatives.
Definition coord.h:90
Analytic equation of state for two fluids (Newtonian case).
double get_kap1() const
Returns the pressure coefficient [unit: ], where .
double get_kap3() const
Returns the pressure coefficient [unit: ], where .
double get_beta() const
Returns the coefficient [unit: ], where .
double get_kap2() const
Returns the pressure coefficient [unit: ], where .
virtual int identify() const =0
Returns a number to identify the sub-classe of Eos_bifluid the object belongs to.
void equilibrium_bi(double ent_c, double ent_c2, double omega0, double omega20, const Tbl &ent_limit, const Tbl &ent2_limit, const Itbl &icontrol, const Tbl &control, Tbl &diff, int mer_mass, double mbar1_wanted, double mbar2_wanted, double aexp_mass)
Computes an equilibrium configuration.
Tenseur gam_euler2
Lorentz factor between the fluid 2 and Eulerian observers
Tenseur enerps_euler
the combination : useful because in the Newtonian limit .
Tenseur sphph_euler
The component of the stress tensor .
virtual double mass_b() const
Total Baryon mass.
Tenseur j_euler
Total angular momentum (flat-space!) 3-vector , which is related to of the "3+1" decomposition,...
virtual void equation_of_state()
Computes the proper baryon and energy densities, as well as pressure and the coefficients Knn,...
double omega2
Rotation angular velocity for fluid 2 ([f_unit] )
const Eos_bifluid & eos
Equation of state for two-fluids model.
double mass_b1() const
Baryon mass of fluid 1.
Tenseur uuu2
Norm of the (fluid no.2) 3-velocity with respect to the eulerian observer.
double mass_b2() const
Baryon mass of fluid 2.
virtual double grv2() const
Error on the virial identity GRV2.
Tenseur ent2
Log-enthalpy for the second fluid.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Definition etoile.h:1625
Tenseur uuu
Norm of u_euler.
Definition etoile.h:1518
double omega
Rotation angular velocity ([f_unit] )
Definition etoile.h:1501
Tenseur & logn
Metric potential = logn_auto.
Definition etoile.h:1521
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition etoile.h:1531
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1560
Tenseur tggg
Metric potential .
Definition etoile.h:1537
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Definition etoile.h:1526
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Definition etoile.h:1608
Tenseur bbb
Metric factor B.
Definition etoile.h:1504
void update_metric()
Computes metric coefficients from known potentials.
Tenseur ak_car
Scalar .
Definition etoile.h:1586
Tenseur & dzeta
Metric potential = beta_auto.
Definition etoile.h:1534
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
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Definition etoile_rot.C:781
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
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Definition etoile.h:1567
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
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1550
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:432
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Definition etoile.h:509
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:471
Map & mp
Mapping associated with the star.
Definition etoile.h:429
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:437
Tenseur shift
Total shift vector.
Definition etoile.h:512
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition etoile.h:468
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition etoile.h:457
Tenseur a_car
Total conformal factor .
Definition etoile.h:515
double ray_pole() const
Coordinate radius at [r_unit].
Basic integer array class.
Definition itbl.h:122
Radial mapping of rather general form.
Definition map.h:2752
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition map_et.C:905
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
Definition map.h:791
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.
Coord sint
Definition map.h:721
virtual void homothetie(double lambda)=0
Sets a new radial scale.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
Definition map.h:807
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.
Coord phi
coordinate centered on the grid
Definition map.h:720
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
Multi-domain grid.
Definition grilles.h:273
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
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:462
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:457
int get_nzone() const
Returns the number of domains.
Definition grilles.h:448
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:452
Parameter storage.
Definition param.h:125
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition param.C:315
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition param.C:1004
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
Definition param.C:453
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition param.C:385
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition param.C:1142
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:246
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition param.C:522
Basic array class.
Definition tbl.h:161
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:361
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition tbl.h:281
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:301
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:824
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tenseur.C:636
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1170
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition tenseur.C:1548
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
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
int get_etat() const
Returns the logical state.
Definition tenseur.h:707
Values and coefficients of a (real-value) function.
Definition valeur.h:287
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition valeur.h:302
void coef() const
Computes the coeffcients of *this.
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:220
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition cmp_math.C:322
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:504
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:348
Cmp cos(const Cmp &)
Cosine.
Definition cmp_math.C:94
Cmp log(const Cmp &)
Neperian logarithm.
Definition cmp_math.C:296
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...
Lorene prototypes.
Definition app_hor.h:64
Standard units of space, time and mass.