LORENE
et_bin_equilibrium.C
1/*
2 * Method of class Etoile_bin to compute an equilibrium configuration
3 *
4 * (see file etoile.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2000-2001 Eric Gourgoulhon
10 * Copyright (c) 2000-2001 Keisuke Taniguchi
11 *
12 * This file is part of LORENE.
13 *
14 * LORENE is free software; you can redistribute it and/or modify
15 * it under the terms of the GNU General Public License as published by
16 * the Free Software Foundation; either version 2 of the License, or
17 * (at your option) any later version.
18 *
19 * LORENE is distributed in the hope that it will be useful,
20 * but WITHOUT ANY WARRANTY; without even the implied warranty of
21 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22 * GNU General Public License for more details.
23 *
24 * You should have received a copy of the GNU General Public License
25 * along with LORENE; if not, write to the Free Software
26 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
27 *
28 */
29
30
31
32
33/*
34 * $Id: et_bin_equilibrium.C,v 1.15 2016/12/05 16:17:52 j_novak Exp $
35 * $Log: et_bin_equilibrium.C,v $
36 * Revision 1.15 2016/12/05 16:17:52 j_novak
37 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
38 *
39 * Revision 1.14 2014/10/13 08:52:55 j_novak
40 * Lorene classes and functions now belong to the namespace Lorene.
41 *
42 * Revision 1.13 2014/10/06 15:13:08 j_novak
43 * Modified #include directives to use c++ syntax.
44 *
45 * Revision 1.12 2009/06/15 09:25:18 k_taniguchi
46 * Improved the rescaling of the domains.
47 *
48 * Revision 1.11 2008/11/14 13:48:06 e_gourgoulhon
49 * Added parameter pent_limit to force the enthalpy values at the
50 * boundaries between the domains, in case of more than one domain inside
51 * the star.
52 *
53 * Revision 1.10 2004/09/28 15:49:23 f_limousin
54 * Improve the rescaling of the domains for nzone = 4 and nzone = 5.
55 *
56 * Revision 1.9 2004/05/13 08:47:01 f_limousin
57 * Decomment the procedure resize.
58 *
59 * Revision 1.8 2004/05/10 10:15:57 f_limousin
60 * Change to avoid a warning in the compilation of Lorene
61 *
62 * Revision 1.7 2004/05/07 12:36:34 f_limousin
63 * Add new member ssjm1_psi in order to have only one function
64 * equilibrium (the same for strange stars and neutron stars)
65 *
66 * Revision 1.6 2004/05/07 08:32:44 k_taniguchi
67 * Introduction of the version without ssjm1_psi.
68 *
69 * Revision 1.5 2004/04/19 11:06:36 f_limousin
70 * Differents call of Etoile_bin::velocity_potential depending on
71 * the equation of state.
72 *
73 * Revision 1.4 2004/03/25 10:29:04 j_novak
74 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
75 *
76 * Revision 1.3 2003/01/17 13:31:13 f_limousin
77 * Add comments
78 *
79 * Revision 1.2 2002/12/10 13:28:03 k_taniguchi
80 * Change the multiplication "*" to "%"
81 * and flat_scalar_prod to flat_scalar_prod_desal.
82 *
83 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
84 * LORENE
85 *
86 * Revision 2.23 2001/08/07 09:43:02 keisuke
87 * Change of the method to set the longest radius of a star
88 * on the first domain.
89 * Addition of a new argument in Etoile_bin::equilibrium : Tbl fact.
90 *
91 * Revision 2.22 2001/06/22 08:54:19 keisuke
92 * Set the inner values of the second domain of ent
93 * by using the outer ones of the first domain.
94 *
95 * Revision 2.21 2001/06/18 12:50:49 keisuke
96 * Addition of the filter for the source of shift vector.
97 *
98 * Revision 2.20 2001/05/17 12:18:47 keisuke
99 * Change of the method to calculate chi from setting position in map
100 * to val_point.
101 *
102 * Revision 2.19 2001/01/16 17:01:53 keisuke
103 * Change the method to set the values on the surface.
104 *
105 * Revision 2.18 2001/01/10 16:42:51 keisuke
106 * Set the inner values of the second domain of logn_auto
107 * by using the outer ones of the first domain.
108 *
109 * Revision 2.17 2000/12/22 13:08:04 eric
110 * precis_adapt = 1e-14 au lieu de 1e-15.
111 * Sorties graphiques commentees.
112 *
113 * Revision 2.16 2000/12/20 10:32:44 eric
114 * Changement important : nz_search = nzet ---> nz_search = nzet + 1
115 *
116 * Revision 2.15 2000/10/23 14:02:16 eric
117 * Modif de Map_et::adapt: on y rentre desormais avec nz_search
118 * (dans le cas present nz_search = nzet).
119 *
120 * Revision 2.14 2000/09/28 12:19:36 keisuke
121 * Construct logn_auto_regu from logn_auto.
122 * This procedure is needed for et_bin_upmetr.C.
123 *
124 * Revision 2.13 2000/05/25 13:48:12 eric
125 * Ajout de l'argument thres_adapt: l'adaptation du mapping n'est
126 * plus effectuee si dH/dr_eq passe sous un certain seuil.
127 *
128 * Revision 2.12 2000/05/25 12:58:31 eric
129 * Modifs classe Param: les int et double sont desormais passes par leurs
130 * adresses.
131 *
132 * Revision 2.11 2000/03/29 11:57:38 eric
133 * *** empty log message ***
134 *
135 * Revision 2.10 2000/03/29 11:53:41 eric
136 * Modif affichage
137 *
138 * Revision 2.9 2000/03/22 16:37:29 eric
139 * Calcul des erreurs dans la resolution des equations de Poisson
140 * et sortie de ces erreurs dans le Tbl diff.
141 *
142 * Revision 2.8 2000/03/22 12:56:18 eric
143 * Nouveau prototype d'Etoile_bin::equilibrium : diff_ent est remplace
144 * par le Tbl diff.
145 *
146 * Revision 2.7 2000/03/10 15:47:19 eric
147 * On appel desormais poisson_vect avec dzpuis = 4.
148 *
149 * Revision 2.6 2000/03/07 16:52:15 eric
150 * Modifs manipulations source pour le shift.
151 *
152 * Revision 2.5 2000/03/07 08:32:47 eric
153 * Appel de Map_radial::reevaluate_sym (pour tenir compte de la symetrie
154 * / plan y=0).
155 *
156 * Revision 2.4 2000/02/17 19:56:57 eric
157 * L'appel de Map_radial::reevaluate pour ent est fait sur nzet+1 zone
158 * et non plus nzet.
159 *
160 * Revision 2.2 2000/02/16 17:12:03 eric
161 * Premiere version avec les equations du champ gravitationnel.
162 *
163 * Revision 2.1 2000/02/15 15:59:52 eric
164 * *** empty log message ***
165 *
166 * Revision 2.0 2000/02/15 15:40:42 eric
167 * *** empty log message ***
168 *
169 *
170 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_equilibrium.C,v 1.15 2016/12/05 16:17:52 j_novak Exp $
171 *
172 */
173
174// Headers C
175#include <cmath>
176
177// Headers Lorene
178#include "etoile.h"
179#include "param.h"
180#include "eos.h"
181
182#include "graphique.h"
183#include "utilitaires.h"
184#include "unites.h"
185
186namespace Lorene {
187void Etoile_bin::equilibrium(double ent_c,
188 int mermax, int mermax_poisson,
189 double relax_poisson, int mermax_potvit,
190 double relax_potvit, double thres_adapt,
191 const Tbl& fact_resize, Tbl& diff, const Tbl* pent_limit ) {
192
193
194 // Fundamental constants and units
195 // -------------------------------
196
197 using namespace Unites ;
198
199 // Initializations
200 // ---------------
201
202 const Mg3d* mg = mp.get_mg() ;
203 int nz = mg->get_nzone() ; // total number of domains
204
205 // The following is required to initialize mp_prev as a Map_et:
206 Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
207
208 // Domain and radial indices of points at the surface of the star:
209 int l_b = nzet - 1 ;
210 int i_b = mg->get_nr(l_b) - 1 ;
211 int k_b ;
212 int j_b ;
213
214 // Value of the enthalpy defining the surface of the star
215 double ent_b = 0 ;
216
217 // Error indicators
218 // ----------------
219
220 double& diff_ent = diff.set(0) ;
221 double& diff_vel_pot = diff.set(1) ;
222 double& diff_logn = diff.set(2) ;
223 double& diff_beta = diff.set(3) ;
224 double& diff_shift_x = diff.set(4) ;
225 double& diff_shift_y = diff.set(5) ;
226 double& diff_shift_z = diff.set(6) ;
227
228 // Parameters for the function Map_et::adapt
229 // -----------------------------------------
230
231 Param par_adapt ;
232 int nitermax = 100 ;
233 int niter ;
234 int adapt_flag = 1 ; // 1 = performs the full computation,
235 // 0 = performs only the rescaling by
236 // the factor alpha_r
237 //## int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
238 int nz_search = nzet ; // Number of domains for searching the enthalpy
239 // isosurfaces
240
241 double precis_secant = 1.e-14 ;
242 double alpha_r ;
243 double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
244
245 par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
246 // locate zeros by the secant method
247 par_adapt.add_int(nzet, 1) ; // number of domains where the adjustment
248 // to the isosurfaces of ent is to be
249 // performed
250 par_adapt.add_int(nz_search, 2) ; // number of domains to search for
251 // the enthalpy isosurface
252 par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
253 // 0 = performs only the rescaling by
254 // the factor alpha_r
255 par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
256 // (theta_*, phi_*)
257 par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
258 // (theta_*, phi_*)
259
260 par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
261 // the secant method
262
263 par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
264 // the determination of zeros by
265 // the secant method
266 par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping, 0 = contracting mapping
267
268 par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
269 // distances will be multiplied
270
271 // Enthalpy values for the adaptation
272 Tbl ent_limit(nzet) ;
273 if (pent_limit != 0x0) ent_limit = *pent_limit ;
274
275 par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
276 // to define the isosurfaces.
277
278 // Parameters for the function Map_et::poisson for logn_auto
279 // ---------------------------------------------------------
280
281 double precis_poisson = 1.e-16 ;
282
283 Param par_poisson1 ;
284
285 par_poisson1.add_int(mermax_poisson, 0) ; // maximum number of iterations
286 par_poisson1.add_double(relax_poisson, 0) ; // relaxation parameter
287 par_poisson1.add_double(precis_poisson, 1) ; // required precision
288 par_poisson1.add_int_mod(niter, 0) ; // number of iterations actually used
289 par_poisson1.add_cmp_mod( ssjm1_logn ) ;
290
291 // Parameters for the function Map_et::poisson for beta_auto
292 // ---------------------------------------------------------
293
294 Param par_poisson2 ;
295
296 par_poisson2.add_int(mermax_poisson, 0) ; // maximum number of iterations
297 par_poisson2.add_double(relax_poisson, 0) ; // relaxation parameter
298 par_poisson2.add_double(precis_poisson, 1) ; // required precision
299 par_poisson2.add_int_mod(niter, 0) ; // number of iterations actually used
300 par_poisson2.add_cmp_mod( ssjm1_beta ) ;
301
302
303 // Parameters for the function Tenseur::poisson_vect
304 // -------------------------------------------------
305
306 Param par_poisson_vect ;
307
308 par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
309 par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
310 par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
311 par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
312 par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
313 par_poisson_vect.add_int_mod(niter, 0) ;
314
315
316 // External potential
317 // See Eq (99) from Gourgoulhon et al. (2001)
318 // -----------------------------------------
319
320 Tenseur pot_ext = logn_comp + pot_centri + loggam ;
321//##
322// des_coupe_z(pot_ext(), 0., 1, "pot_ext", &(ent()) ) ;
323//##
324
325 Tenseur ent_jm1 = ent ; // Enthalpy at previous step
326
327 Tenseur source(mp) ; // source term in the equation for logn_auto
328 // and beta_auto
329
330 Tenseur source_shift(mp, 1, CON, ref_triad) ; // source term in the equation
331 // for shift_auto
332
333 //=========================================================================
334 // Start of iteration
335 //=========================================================================
336
337 for(int mer=0 ; mer<mermax ; mer++ ) {
338
339 cout << "-----------------------------------------------" << endl ;
340 cout << "step: " << mer << endl ;
341 cout << "diff_ent = " << diff_ent << endl ;
342
343 //-----------------------------------------------------
344 // Resolution of the elliptic equation for the velocity
345 // scalar potential
346 //-----------------------------------------------------
347
348 if (irrotational) {
349 diff_vel_pot = velocity_potential(mermax_potvit,
350 precis_poisson, relax_potvit) ;
351
352 }
353
354 //-----------------------------------------------------
355 // Computation of the new radial scale
356 //-----------------------------------------------------
357
358 // alpha_r (r = alpha_r r') is determined so that the enthalpy
359 // takes the requested value ent_b at the stellar surface
360
361 // Values at the center of the star:
362 double logn_auto_c = logn_auto()(0, 0, 0, 0) ;
363 double pot_ext_c = pot_ext()(0, 0, 0, 0) ;
364
365 // Search for the reference point (theta_*, phi_*) [notation of
366 // Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
367 // at the surface of the star
368 // ------------------------------------------------------------
369 double alpha_r2 = 0 ;
370 for (int k=0; k<mg->get_np(l_b); k++) {
371 for (int j=0; j<mg->get_nt(l_b); j++) {
372
373 double pot_ext_b = pot_ext()(l_b, k, j, i_b) ;
374 double logn_auto_b = logn_auto()(l_b, k, j, i_b) ;
375
376
377 // See Eq (100) from Gourgoulhon et al. (2001)
378 double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b) /
379 ( logn_auto_b - logn_auto_c ) ;
380
381// cout << "k, j, alpha_r2_jk : " << k << " " << j << " "
382// << alpha_r2_jk << endl ;
383
384 if (alpha_r2_jk > alpha_r2) {
385 alpha_r2 = alpha_r2_jk ;
386 k_b = k ;
387 j_b = j ;
388 }
389
390 }
391 }
392
393 alpha_r = sqrt(alpha_r2) ;
394
395 cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
396 << alpha_r << endl ;
397
398 // New value of logn_auto
399 // ----------------------
400
401 logn_auto = alpha_r2 * logn_auto ;
402 logn_auto_regu = alpha_r2 * logn_auto_regu ;
403 logn_auto_c = logn_auto()(0, 0, 0, 0) ;
404
405
406 //------------------------------------------------------------
407 // Change the values of the inner points of the second domain
408 // by those of the outer points of the first domain
409 //------------------------------------------------------------
410
411 (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
412
413
414 //------------------------------------------
415 // First integral --> enthalpy in all space
416 // See Eq (98) from Gourgoulhon et al. (2001)
417 //-------------------------------------------
418
419 ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
420
421 (ent().va).smooth(nzet, (ent.set()).va) ;
422
423 //----------------------------------------------------
424 // Adaptation of the mapping to the new enthalpy field
425 //----------------------------------------------------
426
427 // Shall the adaptation be performed (cusp) ?
428 // ------------------------------------------
429
430 double dent_eq = ent().dsdr().val_point(ray_eq(),M_PI/2.,0.) ;
431 double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
432 double rap_dent = fabs( dent_eq / dent_pole ) ;
433 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
434
435 if ( rap_dent < thres_adapt ) {
436 adapt_flag = 0 ; // No adaptation of the mapping
437 cout << "******* FROZEN MAPPING *********" << endl ;
438 }
439 else{
440 adapt_flag = 1 ; // The adaptation of the mapping is to be
441 // performed
442 }
443
444
445 if (pent_limit == 0x0) {
446 ent_limit.set_etat_qcq() ;
447 for (int l=0; l<nzet; l++) { // loop on domains inside the star
448 ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
449 }
450 ent_limit.set(nzet-1) = ent_b ;
451 }
452
453 Map_et mp_prev = mp_et ;
454
455//## cout << "Enthalpy field at the outer boundary of domain 0 : "
456// << endl ;
457// for (int k=0; k<mg->get_np(0); k++) {
458// cout << "k = " << k << " : " ;
459// for (int j=0; j<mg->get_nt(0); j++) {
460// cout << " " << ent()(0, k, j, mg->get_nr(0)-1) ;
461// }
462// cout << endl ;
463// }
464// cout << "Enthalpy field at the inner boundary of domain 1 : "
465// << endl ;
466// for (int k=0; k<mg->get_np(1); k++) {
467// cout << "k = " << k << " : " ;
468// for (int j=0; j<mg->get_nt(1); j++) {
469// cout << " " << ent()(1, k, j, 0) ;
470// }
471// cout << endl ;
472// }
473// cout << "Difference enthalpy field boundary between domains 0 and 1: "
474// << endl ;
475// for (int k=0; k<mg->get_np(1); k++) {
476// cout << "k = " << k << " : " ;
477// for (int j=0; j<mg->get_nt(1); j++) {
478// cout << " " << ent()(0, k, j, mg->get_nr(0)-1) -
479// ent()(1, k, j, 0) ;
480// }
481// cout << endl ;
482// }
483
484
485//##
486// des_coupe_z(gam_euler(), 0., 1, "gam_euler") ;
487// des_coupe_z(loggam(), 0., 1, "loggam") ;
488// des_coupe_y(loggam(), 0., 1, "loggam") ;
489// des_coupe_z(d_psi(0), 0., 1, "d_psi_0") ;
490// des_coupe_z(d_psi(1), 0., 1, "d_psi_1") ;
491// des_coupe_z(d_psi(2), 0., 1, "d_psi_2") ;
492// des_coupe_z(ent(), 0., 1, "ent before adapt", &(ent()) ) ;
493//##
494
495
496 mp.adapt(ent(), par_adapt) ;
497
498 // Readjustment of the external boundary of domain l=nzet
499 // to keep a fixed ratio with respect to star's surface
500
501 double rr_in_1 = mp.val_r(nzet, -1., M_PI/2., 0.) ;
502
503 // Resizes the outer boundary of the shell including the comp. NS
504 double rr_out_nm2 = mp.val_r(nz-2, 1., M_PI/2., 0.) ;
505
506 mp.resize(nz-2, rr_in_1/rr_out_nm2 * fact_resize(1)) ;
507
508 // Resizes the inner boundary of the shell including the comp. NS
509 double rr_out_nm3 = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
510
511 mp.resize(nz-3, rr_in_1/rr_out_nm3 * fact_resize(0)) ;
512
513 if (nz > nzet+3) {
514
515 // Resize of the domain from 1(nzet) to N-4
516 double rr_out_nm3_new = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
517
518 for (int i=nzet-1; i<nz-4; i++) {
519
520 double rr_out_i = mp.val_r(i, 1., M_PI/2., 0.) ;
521
522 double rr_mid = rr_out_i
523 + (rr_out_nm3_new - rr_out_i) / double(nz - 3 - i) ;
524
525 double rr_2timesi = 2. * rr_out_i ;
526
527 if (rr_2timesi < rr_mid) {
528
529 double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
530
531 mp.resize(i+1, rr_2timesi / rr_out_ip1) ;
532
533 }
534 else {
535
536 double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
537
538 mp.resize(i+1, rr_mid / rr_out_ip1) ;
539
540 } // End of else
541
542 } // End of i loop
543
544 } // End of (nz > nzet+3) loop
545
546
547 //----------------------------------------------------
548 // Computation of the enthalpy at the new grid points
549 //----------------------------------------------------
550
551 mp_prev.homothetie(alpha_r) ;
552
553 mp.reevaluate_symy(&mp_prev, nzet+1, ent.set()) ;
554
555// des_coupe_z(ent(), 0., 1, "ent after reevaluate", &(ent()) ) ;
556
557 double ent_s_max = -1 ;
558 int k_s_max = -1 ;
559 int j_s_max = -1 ;
560 for (int k=0; k<mg->get_np(l_b); k++) {
561 for (int j=0; j<mg->get_nt(l_b); j++) {
562 double xx = fabs( ent()(l_b, k, j, i_b) ) ;
563 if (xx > ent_s_max) {
564 ent_s_max = xx ;
565 k_s_max = k ;
566 j_s_max = j ;
567 }
568 }
569 }
570 cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
571 << " and nzet : " << endl ;
572 cout << " " << ent_s_max << " reached for k = " << k_s_max <<
573 " and j = " << j_s_max << endl ;
574
575 //----------------------------------------------------
576 // Equation of state
577 //----------------------------------------------------
578
579 equation_of_state() ; // computes new values for nbar (n), ener (e)
580 // and press (p) from the new ent (H)
581
582 //---------------------------------------------------------
583 // Matter source terms in the gravitational field equations
584 //---------------------------------------------------------
585
586 hydro_euler() ; // computes new values for ener_euler (E),
587 // s_euler (S) and u_euler (U^i)
588
589 //--------------------------------------------------------
590 // Poisson equation for logn_auto (nu_auto)
591 //--------------------------------------------------------
592
593 // Source
594 // See Eq (50) from Gourgoulhon et al. (2001)
595 // ------------------------------------------
596
597 if (relativistic) {
598 source = qpig * a_car % (ener_euler + s_euler)
602
603 }
604 else {
605 source = qpig * nbar ;
606 }
607
608 source.set_std_base() ;
609
610 // Resolution of the Poisson equation
611 // ----------------------------------
612
613 source().poisson(par_poisson1, logn_auto.set()) ;
614
615 // Construct logn_auto_regu for et_bin_upmetr.C
616 // --------------------------------------------
617
619
620 // Check: has the Poisson equation been correctly solved ?
621 // -----------------------------------------------------
622
623 Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
624 cout <<
625 "Relative error in the resolution of the equation for logn_auto : "
626 << endl ;
627 for (int l=0; l<nz; l++) {
628 cout << tdiff_logn(l) << " " ;
629 }
630 cout << endl ;
631 diff_logn = max(abs(tdiff_logn)) ;
632
633
634 if (relativistic) {
635
636 //--------------------------------------------------------
637 // Poisson equation for beta_auto
638 //--------------------------------------------------------
639
640 // Source
641 // See Eq (51) from Gourgoulhon et al. (2001)
642 // ------------------------------------------
643
644 source = qpig * a_car % s_euler
645 + .75 * ( akcar_auto + akcar_comp )
650
651 source.set_std_base() ;
652
653 // Resolution of the Poisson equation
654 // ----------------------------------
655
656 source().poisson(par_poisson2, beta_auto.set()) ;
657
658 // Check: has the Poisson equation been correctly solved ?
659 // -----------------------------------------------------
660
661 Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
662 cout <<
663 "Relative error in the resolution of the equation for beta_auto : "
664 << endl ;
665 for (int l=0; l<nz; l++) {
666 cout << tdiff_beta(l) << " " ;
667 }
668 cout << endl ;
669 diff_beta = max(abs(tdiff_beta)) ;
670
671 //--------------------------------------------------------
672 // Vector Poisson equation for shift_auto
673 //--------------------------------------------------------
674
675 // Source
676 // See Eq (52) from Gourgoulhon et al. (2001)
677 // ------
678
679 Tenseur vtmp = 6. * ( d_beta_auto + d_beta_comp )
680 -8. * ( d_logn_auto + d_logn_comp ) ;
681
682 source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
683 % u_euler
685
686 source_shift.set_std_base() ;
687
688 // Resolution of the Poisson equation
689 // ----------------------------------
690
691 // Filter for the source of shift vector
692
693 for (int i=0; i<3; i++) {
694
695 if (source_shift(i).get_etat() != ETATZERO)
696 source_shift.set(i).filtre(4) ;
697
698 }
699
700 // For Tenseur::poisson_vect, the triad must be the mapping triad,
701 // not the reference one:
702
703 source_shift.change_triad( mp.get_bvect_cart() ) ;
704
705 for (int i=0; i<3; i++) {
706 if(source_shift(i).dz_nonzero()) {
707 assert( source_shift(i).get_dzpuis() == 4 ) ;
708 }
709 else{
710 (source_shift.set(i)).set_dzpuis(4) ;
711 }
712 }
713
714 //##
715 // source_shift.dec2_dzpuis() ; // dzpuis 4 -> 2
716
717 double lambda_shift = double(1) / double(3) ;
718
719 source_shift.poisson_vect(lambda_shift, par_poisson_vect,
721
722
723 // Check: has the equation for shift_auto been correctly solved ?
724 // --------------------------------------------------------------
725
726 // Divergence of shift_auto :
727 Tenseur divn = contract(shift_auto.gradient(), 0, 1) ;
728 divn.dec2_dzpuis() ; // dzpuis 2 -> 0
729
730 // Grad(div) :
731 Tenseur graddivn = divn.gradient() ;
732 graddivn.inc2_dzpuis() ; // dzpuis 2 -> 4
733
734 // Full operator :
735 Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
736 lap_shift.set_etat_qcq() ;
737 for (int i=0; i<3; i++) {
738 lap_shift.set(i) = shift_auto(i).laplacien()
739 + lambda_shift * graddivn(i) ;
740 }
741
742 Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
743 Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
744 Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
745
746 cout <<
747 "Relative error in the resolution of the equation for shift_auto : "
748 << endl ;
749 cout << "x component : " ;
750 for (int l=0; l<nz; l++) {
751 cout << tdiff_shift_x(l) << " " ;
752 }
753 cout << endl ;
754 cout << "y component : " ;
755 for (int l=0; l<nz; l++) {
756 cout << tdiff_shift_y(l) << " " ;
757 }
758 cout << endl ;
759 cout << "z component : " ;
760 for (int l=0; l<nz; l++) {
761 cout << tdiff_shift_z(l) << " " ;
762 }
763 cout << endl ;
764
765 diff_shift_x = max(abs(tdiff_shift_x)) ;
766 diff_shift_y = max(abs(tdiff_shift_y)) ;
767 diff_shift_z = max(abs(tdiff_shift_z)) ;
768
769 // Final result
770 // ------------
771 // The output of Tenseur::poisson_vect is on the mapping triad,
772 // it should therefore be transformed to components on the
773 // reference triad :
774
775 shift_auto.change_triad( ref_triad ) ;
776
777
778 } // End of relativistic equations
779
780
781 if (nzet > 1) {
782 cout.precision(10) ;
783
784 for (int ltrans = 0; ltrans < nzet-1; ltrans++) {
785 cout << endl << "Values at boundary between domains no. " << ltrans << " and " << ltrans+1 << " for theta = pi/2 and phi = 0 :" << endl ;
786
787 double rt1 = mp.val_r(ltrans, 1., M_PI/2, 0.) ;
788 double rt2 = mp.val_r(ltrans+1, -1., M_PI/2, 0.) ;
789 cout << " Coord. r [km] (left, right, rel. diff) : "
790 << rt1 / km << " " << rt2 / km << " " << (rt2 - rt1)/rt1 << endl ;
791
792 int ntm1 = mg->get_nt(ltrans) - 1;
793 int nrm1 = mg->get_nr(ltrans) - 1 ;
794 double ent1 = ent()(ltrans, 0, ntm1, nrm1) ;
795 double ent2 = ent()(ltrans+1, 0, ntm1, 0) ;
796 cout << " Enthalpy (left, right, rel. diff) : "
797 << ent1 << " " << ent2 << " " << (ent2-ent1)/ent1 << endl ;
798
799 double press1 = press()(ltrans, 0, ntm1, nrm1) ;
800 double press2 = press()(ltrans+1, 0, ntm1, 0) ;
801 cout << " Pressure (left, right, rel. diff) : "
802 << press1 << " " << press2 << " " << (press2-press1)/press1 << endl ;
803
804 double nb1 = nbar()(ltrans, 0, ntm1, nrm1) ;
805 double nb2 = nbar()(ltrans+1, 0, ntm1, 0) ;
806 cout << " Baryon density (left, right, rel. diff) : "
807 << nb1 << " " << nb2 << " " << (nb2-nb1)/nb1 << endl ;
808 }
809 }
810/* if (mer % 10 == 0) {
811 cout << "mer = " << mer << endl ;
812 double r_max = 1.2 * ray_eq() ;
813 des_profile(nbar(), 0., r_max, M_PI/2, 0., "n", "Baryon density") ;
814 des_profile(ener(), 0., r_max, M_PI/2, 0., "e", "Energy density") ;
815 des_profile(press(), 0., r_max, M_PI/2, 0., "p", "Pressure") ;
816 des_profile(ent(), 0., r_max, M_PI/2, 0., "H", "Enthalpy") ;
817 }*/
818
819 //-------------------------------------------------
820 // Relative change in enthalpy
821 //-------------------------------------------------
822
823 Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
824 diff_ent = diff_ent_tbl(0) ;
825 for (int l=1; l<nzet; l++) {
826 diff_ent += diff_ent_tbl(l) ;
827 }
828 diff_ent /= nzet ;
829
830
831 ent_jm1 = ent ;
832
833
834 } // End of main loop
835
836 //=========================================================================
837 // End of iteration
838 //=========================================================================
839
840
841}
842}
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition cmp_manip.C:77
Tenseur d_beta_auto
Gradient of beta_auto (Cartesian components with respect to ref_triad ).
Definition etoile.h:882
Tenseur d_logn_auto
Gradient of logn_auto (Cartesian components with respect to ref_triad ).
Definition etoile.h:862
Tenseur akcar_comp
Part of the scalar generated by shift_auto and shift_comp , i.e.
Definition etoile.h:947
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:976
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
void equilibrium(double ent_c, int mermax, int mermax_poisson, double relax_poisson, int mermax_potvit, double relax_potvit, double thres_adapt, const Tbl &fact, Tbl &diff, const Tbl *ent_limit=0x0)
Computes an equilibrium configuration.
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur 's...
Definition etoile.h:831
bool irrotational
true for an irrotational star, false for a corotating one
Definition etoile.h:825
Tenseur w_shift
Vector used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:911
Tenseur logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principaly by ...
Definition etoile.h:857
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad ).
Definition etoile.h:872
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 by means...
Definition etoile.h:986
Cmp ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto by mea...
Definition etoile.h:962
Tenseur pot_centri
Centrifugal potential.
Definition etoile.h:956
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition etoile.h:852
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition etoile.h:892
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor generated by shift_auto .
Definition etoile.h:928
Tenseur akcar_auto
Part of the scalar generated by shift_auto , i.e.
Definition etoile.h:941
Cmp ssjm1_beta
Effective source at the previous step for the resolution of the Poisson equation for beta_auto by mea...
Definition etoile.h:968
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad ).
Definition etoile.h:887
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:921
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:435
Tenseur logn_auto_regu
Regular part of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:494
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Definition etoile.h:512
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:487
Tenseur nbar
Baryon density in the fluid frame.
Definition etoile.h:462
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition etoile.C:569
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition etoile.h:477
Map & mp
Mapping associated with the star.
Definition etoile.h:432
Tenseur press
Fluid pressure.
Definition etoile.h:464
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:440
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition etoile.h:468
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition etoile.h:471
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Definition etoile.h:460
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition etoile.h:509
Tenseur a_car
Total conformal factor .
Definition etoile.h:518
double ray_pole() const
Coordinate radius at [r_unit].
Radial mapping of rather general form.
Definition map.h:2770
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition map_et.C:928
Multi-domain grid.
Definition grilles.h:279
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:479
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:474
int get_nzone() const
Returns the number of domains.
Definition grilles.h:465
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:469
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:318
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition param.C:1007
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition param.C:388
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition param.C:1145
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:249
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition param.C:525
Basic array class.
Definition tbl.h:161
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:364
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:304
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:830
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tenseur.C:642
void dec2_dzpuis()
dzpuis -= 2 ;
Definition tenseur.C:1136
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates).
Definition tenseur.C:1548
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1176
void inc2_dzpuis()
dzpuis += 2 ;
Definition tenseur.C:1149
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition tenseur.C:674
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:507
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition cmp_math.C:438
Cmp abs(const Cmp &)
Absolute value.
Definition cmp_math.C:413
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tenseur flat_scalar_prod_desal(const Tenseur &t1, const Tenseur &t2)
Same as flat_scalar_prod but with desaliasing.
Lorene prototypes.
Definition app_hor.h:67
virtual void laplacien(const Scalar &uu, int zec_mult_r, Scalar &lap) const =0
Computes the Laplacian of a scalar field.
Standard units of space, time and mass.