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