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