LORENE
star_bhns_equilibrium.C
1/*
2 * Method of class Star_bhns to compute an equilibrium configuration
3 * of a neutron star in a black hole-neutron star binary
4 *
5 * (see file star_bhns.h for documentation).
6 *
7 */
8
9/*
10 * Copyright (c) 2005-2007 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 version 2
16 * as published by the Free Software Foundation.
17 *
18 * LORENE is distributed in the hope that it will be useful,
19 * but WITHOUT ANY WARRANTY; without even the implied warranty of
20 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 * GNU General Public License for more details.
22 *
23 * You should have received a copy of the GNU General Public License
24 * along with LORENE; if not, write to the Free Software
25 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26 *
27 */
28
29
30
31/*
32 * $Id: star_bhns_equilibrium.C,v 1.6 2018/11/16 14:34:37 j_novak Exp $
33 * $Log: star_bhns_equilibrium.C,v $
34 * Revision 1.6 2018/11/16 14:34:37 j_novak
35 * Changed minor points to avoid some compilation warnings.
36 *
37 * Revision 1.5 2016/12/05 16:18:16 j_novak
38 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
39 *
40 * Revision 1.4 2014/10/13 08:53:40 j_novak
41 * Lorene classes and functions now belong to the namespace Lorene.
42 *
43 * Revision 1.3 2014/10/06 15:13:16 j_novak
44 * Modified #include directives to use c++ syntax.
45 *
46 * Revision 1.2 2008/05/15 19:13:45 k_taniguchi
47 * Change of some parameters.
48 *
49 * Revision 1.1 2007/06/22 01:30:45 k_taniguchi
50 * *** empty log message ***
51 *
52 *
53 * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_equilibrium.C,v 1.6 2018/11/16 14:34:37 j_novak Exp $
54 *
55 */
56
57// C++ headers
58//#include <>
59
60// C headers
61#include <cmath>
62
63// Lorene headers
64#include "star_bhns.h"
65#include "param.h"
66#include "cmp.h"
67#include "tenseur.h"
68#include "utilitaires.h"
69#include "unites.h"
70//#include "graphique.h"
71
72namespace Lorene {
73void Star_bhns::equilibrium_bhns(double ent_c, const double& mass_bh,
74 const double& sepa, bool kerrschild,
75 int, int mermax_ns, int mermax_potvit,
76 int mermax_poisson, int filter_r,
77 int filter_r_s, int filter_p_s,
78 double relax_poisson, double relax_potvit,
79 double thres_adapt, double resize_ns,
80 const Tbl& fact_resize, Tbl& diff) {
81
82 // Fundamental constants and units
83 // -------------------------------
84 using namespace Unites ;
85
86 // Initializations
87 // ---------------
88
89 const Mg3d* mg = mp.get_mg() ;
90 int nz = mg->get_nzone() ; // Total number of domain
91 // int nt = mg->get_nt(0) ;
92
93 // The following is required to initialize mp_prev as a Map_et:
94 Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
95
96 // Domain and radial indices of points at the surface of the star:
97 int l_b = nzet - 1 ;
98 int i_b = mg->get_nr(l_b) - 1 ;
99 int k_b ;
100 int j_b ;
101
102 // Value of the enthalpy defining the surface of the star
103 double ent_b = 0 ;
104
105 // Error indicators
106 // ----------------
107
108 double& diff_ent = diff.set(0) ;
109 double& diff_vel_pot = diff.set(1) ;
110 double& diff_lapconf = diff.set(2) ;
111 double& diff_confo = diff.set(3) ;
112 double& diff_shift_x = diff.set(4) ;
113 double& diff_shift_y = diff.set(5) ;
114 double& diff_shift_z = diff.set(6) ;
115 double& diff_dHdr = diff.set(7) ;
116 double& diff_dHdr_min = diff.set(8) ;
117 double& diff_phi_min = diff.set(9) ;
118 double& diff_radius = diff.set(10) ;
119
120 // Parameters for the function Map_et::adapt
121 // -----------------------------------------
122
123 Param par_adapt ;
124 int nitermax = 100 ;
125 int niter ;
126 int adapt_flag = 1 ; // 1 = performs the full computation,
127 // 0 = performs only the rescaling by
128 // the factor alpha_r
129
130 // Number of domains for searching the enthalpy isosurfaces
131 int nz_search = nzet + 1 ;
132 // int nz_search = nzet ;
133
134 double precis_secant = 1.e-14 ;
135 double alpha_r ;
136 double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
137
138 Tbl ent_limit(nz) ;
139
140 par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
141 // locate zeros by the secant method
142 par_adapt.add_int(nzet, 1) ; // number of domains where the
143 // adjustment to the isosurfaces of
144 // ent is to be performed
145 par_adapt.add_int(nz_search, 2) ; // number of domains to search for
146 // the enthalpy isosurface
147 par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
148 // 0 = performs only the rescaling by
149 // the factor alpha_r
150 par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
151 // (theta_*, phi_*)
152 par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
153 // (theta_*, phi_*)
154 par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used
155 // in the secant method
156 par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
157 // the determination of zeros by
158 // the secant method
159 par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping,
160 // 0 = contracting mapping
161 par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
162 // distances will be multiplied
163 par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
164 // to define the isosurfaces
165
166 Cmp ssjm1lapconf(ssjm1_lapconf) ;
167 Cmp ssjm1confo(ssjm1_confo) ;
168
169 ssjm1lapconf.set_etat_qcq() ;
170 ssjm1confo.set_etat_qcq() ;
171
172 double precis_poisson = 1.e-14 ;
173
174 // Parameters for the function Scalar::poisson for lapconf_auto
175 // ------------------------------------------------------------
176
177 Param par_lapconf ;
178
179 par_lapconf.add_int(mermax_poisson, 0) ; // maximum number of iterations
180 par_lapconf.add_double(relax_poisson, 0) ; // relaxation parameter
181 par_lapconf.add_double(precis_poisson, 1) ; // required precision
182 par_lapconf.add_int_mod(niter, 0) ; // number of iterations actually used
183 par_lapconf.add_cmp_mod( ssjm1lapconf ) ;
184
185 // Parameters for the function Scalar::poisson for confo_auto
186 // ----------------------------------------------------------
187
188 Param par_confo ;
189
190 par_confo.add_int(mermax_poisson, 0) ; // maximum number of iterations
191 par_confo.add_double(relax_poisson, 0) ; // relaxation parameter
192 par_confo.add_double(precis_poisson, 1) ; // required precision
193 par_confo.add_int_mod(niter, 0) ; // number of iterations actually used
194 par_confo.add_cmp_mod( ssjm1confo ) ;
195
196 // Parameters for the function Vector::poisson for shift method 2
197 // --------------------------------------------------------------
198
199 Param par_shift2 ;
200
201 par_shift2.add_int(mermax_poisson, 0) ; // maximum number of iterations
202 par_shift2.add_double(relax_poisson, 0) ; // relaxation parameter
203 par_shift2.add_double(precis_poisson, 1) ; // required precision
204 par_shift2.add_int_mod(niter, 0) ; // number of iterations actually used
205
206 Cmp ssjm1khi(ssjm1_khi) ;
207 ssjm1khi.set_etat_qcq() ;
208
209 Tenseur ssjm1wshift(mp, 1, CON, mp.get_bvect_cart()) ;
210 ssjm1wshift.set_etat_qcq() ;
211 for (int i=0; i<3; i++) {
212 ssjm1wshift.set(i) = Cmp(ssjm1_wshift(i+1)) ;
213 }
214
215 par_shift2.add_cmp_mod(ssjm1khi) ;
216 par_shift2.add_tenseur_mod(ssjm1wshift) ;
217
218 Scalar ent_jm1 = ent ; // Enthalpy at previous step
219
220 Scalar lapconf_m1(mp) ; // = lapconf_auto - 0.5
221 Scalar confo_m1(mp) ; // = confo_auto - 0.5
222
223 Scalar source_lapconf(mp) ; // Source term in the equation for lapconf_auto
224 source_lapconf.set_etat_qcq() ;
225 Scalar source_confo(mp) ; // Source term in the equation for confo_auto
226 source_confo.set_etat_qcq() ;
227 Vector source_shift(mp, CON, mp.get_bvect_cart()) ; // Source term
228 // in the equation for shift_auto
229 source_shift.set_etat_qcq() ;
230
231 //===========================================================
232 // Start of iteration
233 //===========================================================
234
235 for (int mer_ns=0; mer_ns<mermax_ns; mer_ns++) {
236
237 cout << "-----------------------------------------------" << endl ;
238 cout << "step: " << mer_ns << endl ;
239 cout << "diff_ent = " << diff_ent << endl ;
240
241 //------------------------------------------------------
242 // Resolution of the elliptic equation for the velocity
243 // scalar potential
244 //------------------------------------------------------
245
246 if (irrotational) {
247 diff_vel_pot = velo_pot_bhns(mass_bh, sepa, kerrschild,
248 mermax_potvit, precis_poisson,
249 relax_potvit) ;
250 }
251 else {
252 diff_vel_pot = 0. ; // to avoid the warning
253 }
254
255 //-------------------------------------
256 // Computation of the new radial scale
257 //-------------------------------------
258
259 // alpha_r (r = alpha_r r') is determined so that the enthalpy
260 // takes the requested value ent_b at the stellar surface
261
262 // Values at the center of the star:
263 // lapconf_auto = lapconf_auto=0(r->infty) + 0.5
264 // The term lapconf_auto=0(r->infty) should be rescaled.
265 double lapconf_auto_c = lapconf_auto.val_grid_point(0,0,0,0) - 0.5 ;
266 double lapconf_comp_c = lapconf_comp.val_grid_point(0,0,0,0) ;
267
268 double confo_c = confo_tot.val_grid_point(0,0,0,0) ;
269
270 double gam_c = gam.val_grid_point(0,0,0,0) ;
271 double gam0_c = gam0.val_grid_point(0,0,0,0) ;
272
273 double hhh_c = exp(ent_c) ;
274 double hhh_b = exp(ent_b) ;
275
276 // Search for the reference point (theta_*, phi_*) [notation of
277 // Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
278 // at the surface of the star
279 // ------------------------------------------------------------
280 double alpha_r2 = 0. ;
281
282 for (int k=0; k<mg->get_np(l_b); k++) {
283 for (int j=0; j<mg->get_nt(l_b); j++) {
284
285 double lapconf_auto_b =
286 lapconf_auto.val_grid_point(l_b,k,j,i_b) - 0.5 ;
287 double lapconf_comp_b =
288 lapconf_comp.val_grid_point(l_b,k,j,i_b) ;
289
290 double confo_b = confo_tot.val_grid_point(l_b,k,j,i_b) ;
291
292 double gam_b = gam.val_grid_point(l_b,k,j,i_b) ;
293 double gam0_b = gam0.val_grid_point(l_b,k,j,i_b) ;
294
295 double aaa = (gam0_c*gam_b*hhh_b*confo_c)
296 / (gam0_b*gam_c*hhh_c*confo_b) ;
297
298 // See Eq (100) from Gourgoulhon et al. (2001)
299 double alpha_r2_jk = (aaa * lapconf_comp_b - lapconf_comp_c
300 + 0.5 * (aaa - 1.))
301 / (lapconf_auto_c - aaa * lapconf_auto_b ) ;
302
303 if (alpha_r2_jk > alpha_r2) {
304 alpha_r2 = alpha_r2_jk ;
305 k_b = k ;
306 j_b = j ;
307 }
308 }
309 }
310
311 alpha_r = sqrt(alpha_r2) ;
312
313 cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
314 << alpha_r << endl ;
315
316 // New value of lapconf_auto
317 // -------------------------
318
319 lapconf_auto = alpha_r2 * (lapconf_auto - 0.5) + 0.5 ;
320 Scalar lapconf_tot_tmp = lapconf_auto + lapconf_comp ;
321 lapconf_tot_tmp.std_spectral_base() ;
322
323 /*
324 confo_auto = alpha_r2 * (confo_auto - 0.5) + 0.5 ;
325 Scalar confo_tot_tmp = confo_auto + confo_comp ;
326 confo_tot_tmp.std_spectral_base() ;
327 */
328 //------------------------------------------------------------
329 // Change the values of the inner points of the second domain
330 // by those of the outer points of the first domain
331 //------------------------------------------------------------
332
333 lapconf_auto.set_spectral_va().smooth(nzet,lapconf_auto.set_spectral_va()) ;
334
335 //--------------------------------------------
336 // First integral --> enthalpy in all space
337 // See Eq (98) from Gourgoulhon et al. (2001)
338 //--------------------------------------------
339
340 double log_lapconf_c = log(lapconf_tot_tmp.val_grid_point(0,0,0,0)) ;
341 double log_confo_c = log(confo_tot.val_grid_point(0,0,0,0)) ;
342 double loggam_c = loggam.val_grid_point(0,0,0,0) ;
343 double pot_centri_c = pot_centri.val_grid_point(0,0,0,0) ;
344
345 ent = (ent_c + log_lapconf_c - log_confo_c + loggam_c + pot_centri_c)
346 - log(lapconf_tot_tmp) + log(confo_tot) - loggam - pot_centri ;
347 ent.std_spectral_base() ;
348
349
350 //----------------------------------------------------------
351 // Change the enthalpy field to be set its maximum position
352 // at the coordinate center
353 //----------------------------------------------------------
354
355 double dentdx = ent.dsdx().val_grid_point(0,0,0,0) ;
356 double dentdy = ent.dsdy().val_grid_point(0,0,0,0) ;
357
358 cout << "dH/dx|_center = " << dentdx << endl ;
359 cout << "dH/dy|_center = " << dentdy << endl ;
360
361 double dec_fact_x = 1. ;
362 double dec_fact_y = 1. ;
363
364 Scalar func_in(mp) ;
365 func_in = 1. - dec_fact_x * (dentdx/ent_c) * mp.x
366 - dec_fact_y * (dentdy/ent_c) * mp.y ;
367
368 func_in.annule(nzet, nz-1) ;
369 func_in.std_spectral_base() ;
370
371 Scalar func_ex(mp) ;
372 func_ex = 1. ;
373 func_ex.annule(0, nzet-1) ;
374 func_ex.std_spectral_base() ;
375
376 // New enthalpy field
377 // ------------------
378 ent = ent * (func_in + func_ex) ;
379
380 (ent.set_spectral_va()).smooth(nzet, ent.set_spectral_va()) ;
381
382 double dentdx_new = ent.dsdx().val_grid_point(0,0,0,0) ;
383 double dentdy_new = ent.dsdy().val_grid_point(0,0,0,0) ;
384 cout << "dH/dx|_new = " << dentdx_new << endl ;
385 cout << "dH/dy|_new = " << dentdy_new << endl ;
386
387 //-----------------------------------------------------
388 // Adaptation of the mapping to the new enthalpy field
389 //-----------------------------------------------------
390
391 double dent_eq = ent.dsdr().val_point(ray_eq_pi(),M_PI/2.,M_PI) ;
392 double dent_pole = ent.dsdr().val_point(ray_pole(),0.,0.) ;
393 double rap_dent = fabs( dent_eq / dent_pole ) ;
394 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
395 diff_dHdr = rap_dent ;
396
397 if ( rap_dent < thres_adapt ) {
398 adapt_flag = 0 ; // No adaptation of the mapping
399 cout << "******* FROZEN MAPPING *********" << endl ;
400 }
401 else {
402 adapt_flag = 1 ; // The adaptation of the mapping is to be
403 // performed
404 }
405
406 ent_limit.set_etat_qcq() ;
407 for (int l=0; l<nzet; l++) { // loop on domains inside the star
408 ent_limit.set(l) = ent.val_grid_point(l, k_b, j_b, i_b) ;
409 }
410 ent_limit.set(nzet-1) = ent_b ;
411
412 Map_et mp_prev = mp_et ;
413
414 Cmp ent_cmp(ent) ;
415 mp.adapt(ent_cmp, par_adapt) ;
416 ent = ent_cmp ;
417
418 // Readjustment of the external boundary of domain l=nzet
419 // to keep a fixed ratio with respect to star's surface
420
421 double rr_in_1 = mp.val_r(1, -1., M_PI/2., 0.) ;
422
423 // Resize of the outer boundary of the shell including the BH
424 double rr_out_nm2 = mp.val_r(nz-2, 1., M_PI/2., 0.) ;
425 mp.resize(nz-2, rr_in_1/rr_out_nm2 * fact_resize(1)) ;
426
427 // Resize of the inner boundary of the shell including the BH
428 double rr_out_nm3 = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
429 mp.resize(nz-3, rr_in_1/rr_out_nm3 * fact_resize(0)) ;
430
431 // if (mer % 2 == 0) {
432
433 if (nz > 4) {
434
435 // Resize of the domain 1
436 double rr_out_1 = mp.val_r(1, 1., M_PI/2., 0.) ;
437 mp.resize(1, rr_in_1/rr_out_1 * resize_ns) ;
438
439 if (nz > 5) {
440
441 // Resize of the domain from 2 to N-4
442 double rr_out_nm3_new = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
443
444 for (int i=1; i<nz-4; i++) {
445
446 double rr_out_i = mp.val_r(i, 1., M_PI/2., 0.) ;
447
448 double rr_mid = rr_out_i
449 + (rr_out_nm3_new - rr_out_i) / double(nz - 3 - i) ;
450
451 double rr_2timesi = 2. * rr_out_i ;
452
453 if (rr_2timesi < rr_mid) {
454
455 double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
456 mp.resize(i+1, rr_2timesi / rr_out_ip1) ;
457
458 }
459 else {
460
461 double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
462 mp.resize(i+1, rr_mid / rr_out_ip1) ;
463
464 } // End of else
465
466 } // End of i loop
467
468 } // End of (nz > 5) loop
469
470 } // End of (nz > 4) loop
471
472 // } // End of (mer % 2) loop
473
474 //----------------------------------------------------
475 // Computation of the enthalpy at the new grid points
476 //----------------------------------------------------
477
478 mp_prev.homothetie(alpha_r) ;
479
480 Cmp ent_cmp2 (ent) ;
481 mp.reevaluate(&mp_prev, nzet+1, ent_cmp2) ;
482 ent = ent_cmp2 ;
483
484 double ent_s_max = -1 ;
485 int k_s_max = -1 ;
486 int j_s_max = -1 ;
487
488 for (int k=0; k<mg->get_np(l_b); k++) {
489 for (int j=0; j<mg->get_nt(l_b); j++) {
490 double xx = fabs( ent.val_grid_point(l_b, k, j, i_b) ) ;
491 if (xx > ent_s_max) {
492 ent_s_max = xx ;
493 k_s_max = k ;
494 j_s_max = j ;
495 }
496 }
497 }
498 cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
499 << " and nzet : " << endl ;
500 cout << " " << ent_s_max << " reached for k = " << k_s_max
501 << " and j = " << j_s_max << endl ;
502
503 //-------------------
504 // Equation of state
505 //-------------------
506
507 equation_of_state() ; // computes new values for nbar (n), ener (e)
508 // and press (p) from the new ent (H)
509
510 //----------------------------------------------------------
511 // Matter source terms in the gravitational field equations
512 //----------------------------------------------------------
513
514 hydro_euler_bhns(kerrschild, mass_bh, sepa) ;
515 // computes new values for ener_euler (E),
516 // s_euler (S) and u_euler (U^i)
517
518
519 //-------------------------------------------------
520 // Computation of the minimum of the indicator chi
521 //-------------------------------------------------
522 double azimu_min = phi_min() ;
523 double rad_chi_min = radius_p(azimu_min) ;
524 double chi_min = chi_rp(rad_chi_min, azimu_min) ;
525
526 cout << "| dH/dr_eq / dH/dr_pole | (minimum) = " << chi_min << endl ;
527 cout << " phi = " << azimu_min / M_PI << " [M_PI]" << endl ;
528 cout << " radius = " << rad_chi_min / km << " [km]" << endl ;
529
530 diff_dHdr_min = chi_min ;
531 diff_phi_min = azimu_min ;
532 diff_radius = rad_chi_min ;
533
534 //-----------------------------------
535 // Poisson equation for lapconf_auto
536 //-----------------------------------
537
538 // Source
539 //--------
540
541 Scalar sou_lap1(mp) ; // dzpuis = 0
542 sou_lap1 = qpig * lapconf_tot_tmp * pow(confo_tot,4.)
543 * (0.5*ener_euler + s_euler) ;
544
545 sou_lap1.std_spectral_base() ;
546 sou_lap1.annule(nzet,nz-1) ;
547 sou_lap1.inc_dzpuis(4) ; // dzpuis : 0 -> 4
548
549 Scalar sou_lap2(mp) ; // dzpuis = 4
550 sou_lap2 = 0.875 * (lapconf_auto+0.5) * taij_quad_auto
551 / pow(confo_auto+0.5,8.) ;
552 sou_lap2.std_spectral_base() ;
553
554 source_lapconf = sou_lap1 + sou_lap2 ;
555
556 source_lapconf.std_spectral_base() ;
557 // source_lapse.annule(nzet,nz-1) ;
558
559 if (filter_r != 0) {
560 if (source_lapconf.get_etat() != ETATZERO) {
561 source_lapconf.filtre(filter_r) ;
562 // source_lapse.filtre_r(filter_r,0) ;
563 }
564 }
565
566 assert(source_lapconf.get_dzpuis() == 4) ;
567
568 // Resolution of the Poisson equation (Outer BC : lapconf_m1 -> 0)
569 // ----------------------------------
570
571 lapconf_m1.set_etat_qcq() ;
572 lapconf_m1 = lapconf_auto - 0.5 ;
573 source_lapconf.poisson(par_lapconf, lapconf_m1) ;
574 ssjm1_lapconf = ssjm1lapconf ;
575
576 // Check: has the Poisson equation been correctly solved ?
577 // -------------------------------------------------------
578
579 Tbl tdiff_lapconf = diffrel(lapconf_m1.laplacian(), source_lapconf) ;
580 cout <<
581 "Relative error in the resolution of the equation for lapconf_auto : "
582 << endl ;
583 for (int l=0; l<nz; l++) {
584 cout << tdiff_lapconf(l) << " " ;
585 }
586 cout << endl ;
587 diff_lapconf = max(abs(tdiff_lapconf)) ;
588
589 // Re-construction of the lapconf function
590 // ---------------------------------------
591 lapconf_auto = lapconf_m1 + 0.5 ;
592 // lapconf_tot = lapconf_auto + lapconf_comp
593 // lapconf_auto, _comp -> 0.5 (r -> inf)
594 // lapconf_tot -> 1 (r -> inf)
595
596 //---------------------------------
597 // Poisson equation for confo_auto
598 //---------------------------------
599
600 // Source
601 //--------
602
603 Scalar sou_con1(mp) ; // dzpuis = 0
604 sou_con1 = - 0.5 * qpig * pow(confo_tot,5.) * ener_euler ;
605 sou_con1.std_spectral_base() ;
606 sou_con1.annule(nzet,nz-1) ;
607 sou_con1.inc_dzpuis(4) ; // dzpuis : 0 -> 4
608
609 Scalar sou_con2(mp) ; // dzpuis = 4
610 sou_con2 = - 0.125 * taij_quad_auto / pow(confo_auto+0.5,7.) ;
611 sou_con2.std_spectral_base() ;
612
613 source_confo = sou_con1 + sou_con2 ;
614
615 source_confo.std_spectral_base() ;
616 // source_confo.annule(nzet,nz-1) ;
617
618 if (filter_r != 0) {
619 if (source_confo.get_etat() != ETATZERO) {
620 source_confo.filtre(filter_r) ;
621 // source_confo.filtre_r(filter_r,0) ;
622 }
623 }
624
625 assert(source_confo.get_dzpuis() == 4) ;
626
627 // Resolution of the Poisson equation (Outer BC : confo_m1 -> 0)
628 // ----------------------------------
629
630 confo_m1.set_etat_qcq() ;
631 confo_m1 = confo_auto - 0.5 ;
632 source_confo.poisson(par_confo, confo_m1) ;
633 ssjm1_confo = ssjm1confo ;
634
635 // Check: has the Poisson equation been correctly solved ?
636 // -------------------------------------------------------
637
638 Tbl tdiff_confo = diffrel(confo_m1.laplacian(), source_confo) ;
639 cout <<
640 "Relative error in the resolution of the equation for confo_auto : "
641 << endl ;
642 for (int l=0; l<nz; l++) {
643 cout << tdiff_confo(l) << " " ;
644 }
645 cout << endl ;
646 diff_confo = max(abs(tdiff_confo)) ;
647
648 // Re-construction of the conformal factor
649 // ---------------------------------------
650 confo_auto = confo_m1 + 0.5 ; // confo_tot = confo_auto + confo_comp
651 // confo_auto, _comp -> 0.5 (r -> inf)
652 // confo_tot -> 1 (r -> inf)
653
654 //----------------------------------------
655 // Vector Poisson equation for shift_auto
656 //----------------------------------------
657
658 // Source
659 // ------
660
661 Vector sou_shif1(mp, CON, mp.get_bvect_cart()) ; // dzpuis = 0
662 sou_shif1.set_etat_qcq() ;
663
664 for (int i=1; i<=3; i++) {
665 sou_shif1.set(i) = 4.*qpig * lapconf_tot_tmp
666 * pow(confo_tot, 3.)
667 * (ener_euler + press) * u_euler(i) ;
668 }
669
670 sou_shif1.std_spectral_base() ;
671 sou_shif1.annule(nzet, nz-1) ;
672
673 for (int i=1; i<=3; i++) {
674 (sou_shif1.set(i)).inc_dzpuis(4) ; // dzpuis: 0 -> 4
675 }
676
677 Vector sou_shif2(mp, CON, mp.get_bvect_cart()) ; // dzpuis = 4
678 sou_shif2.set_etat_qcq() ;
679 for (int i=1; i<=3; i++) {
680 sou_shif2.set(i) = 2. *
681 (taij_auto(i,1)*(d_lapconf_auto(1)
682 -7.*(lapconf_auto+0.5)*d_confo_auto(1)
683 /(confo_auto+0.5))
684 +taij_auto(i,2)*(d_lapconf_auto(2)
685 -7.*(lapconf_auto+0.5)*d_confo_auto(2)
686 /(confo_auto+0.5))
687 +taij_auto(i,3)*(d_lapconf_auto(3)
688 -7.*(lapconf_auto+0.5)*d_confo_auto(3)
689 /(confo_auto+0.5))
690 ) / pow(confo_auto+0.5,7.) ;
691 }
692 sou_shif2.std_spectral_base() ;
693
694 source_shift = sou_shif1 + sou_shif2 ;
695
696 source_shift.std_spectral_base() ;
697 // source_shift.annule(nzet, nz-1) ;
698
699 // Resolution of the Poisson equation
700 // ----------------------------------
701
702 // Filter for the source of shift vector
703
704 if (filter_r_s != 0) {
705 for (int i=1; i<=3; i++) {
706 if (source_shift(i).get_etat() != ETATZERO) {
707 source_shift.set(i).filtre(filter_r_s) ;
708 // source_shift.set(i).filtre_r(filter_r_s, 0) ;
709 }
710 }
711 }
712
713 if (filter_p_s != 0) {
714 for (int i=1; i<=3; i++) {
715 if (source_shift(i).get_etat() != ETATZERO) {
716 (source_shift.set(i)).filtre_phi(filter_p_s, nz-1) ;
717 // (source_shift.set(i)).filtre_phi(filter_p_s, 0) ;
718 }
719 }
720 }
721
722 for (int i=1; i<=3; i++) {
723 if(source_shift(i).dz_nonzero()) {
724 assert( source_shift(i).get_dzpuis() == 4 ) ;
725 }
726 else {
727 (source_shift.set(i)).set_dzpuis(4) ;
728 }
729 }
730
731 double lambda = double(1) / double(3) ;
732
733 Tenseur source_p(mp, 1, CON, mp.get_bvect_cart() ) ;
734 source_p.set_etat_qcq() ;
735 for (int i=0; i<3; i++) {
736 source_p.set(i) = Cmp(source_shift(i+1)) ;
737 }
738
739 Tenseur vect_auxi(mp, 1, CON, mp.get_bvect_cart()) ;
740 vect_auxi.set_etat_qcq() ;
741 for (int i=0; i<3 ;i++) {
742 vect_auxi.set(i) = 0. ;
743 }
744 Tenseur scal_auxi(mp) ;
745 scal_auxi.set_etat_qcq() ;
746 scal_auxi.set().annule_hard() ;
747 scal_auxi.set_std_base() ;
748
749 Tenseur resu_p(mp, 1, CON, mp.get_bvect_cart() ) ;
750 resu_p.set_etat_qcq() ;
751
752 source_p.poisson_vect(lambda, par_shift2, resu_p, vect_auxi,
753 scal_auxi) ;
754
755 for (int i=1; i<=3; i++)
756 shift_auto.set(i) = resu_p(i-1) ;
757
758 ssjm1_khi = ssjm1khi ;
759
760 for (int i=0; i<3; i++) {
761 ssjm1_wshift.set(i+1) = ssjm1wshift(i) ;
762 }
763
764 // Check: has the equation for shift_auto been correctly solved ?
765 // --------------------------------------------------------------
766
767 Vector lap_shift = shift_auto.derive_con(flat).divergence(flat)
768 + lambda * shift_auto.divergence(flat).derive_con(flat) ;
769
770 source_shift.dec_dzpuis() ;
771
772 Tbl tdiff_shift_x = diffrel(lap_shift(1), source_shift(1)) ;
773 Tbl tdiff_shift_y = diffrel(lap_shift(2), source_shift(2)) ;
774 Tbl tdiff_shift_z = diffrel(lap_shift(3), source_shift(3)) ;
775
776 cout <<
777 "Relative error in the resolution of the equation for shift_auto : "
778 << endl ;
779 cout << "x component : " ;
780 for (int l=0; l<nz; l++) {
781 cout << tdiff_shift_x(l) << " " ;
782 }
783 cout << endl ;
784 cout << "y component : " ;
785 for (int l=0; l<nz; l++) {
786 cout << tdiff_shift_y(l) << " " ;
787 }
788 cout << endl ;
789 cout << "z component : " ;
790 for (int l=0; l<nz; l++) {
791 cout << tdiff_shift_z(l) << " " ;
792 }
793 cout << endl ;
794
795 diff_shift_x = max(abs(tdiff_shift_x)) ;
796 diff_shift_y = max(abs(tdiff_shift_y)) ;
797 diff_shift_z = max(abs(tdiff_shift_z)) ;
798
799
800 //-----------------------------
801 // Relative change in enthalpy
802 //-----------------------------
803
804 Tbl diff_ent_tbl = diffrel( ent, ent_jm1 ) ;
805 diff_ent = diff_ent_tbl(0) ;
806 for (int l=0; l<nzet; l++) {
807 diff_ent += diff_ent_tbl(l) ;
808 }
809 diff_ent /= nzet ;
810
811 ent_jm1 = ent ;
812
813 /*
814 des_profile( lapconf_auto, 0., 10.,
815 M_PI/2., M_PI, "Self lapconf function of NS",
816 "Lapconf (theta=pi/2, phi=0)" ) ;
817
818 des_profile( lapconf_tot, 0., 10.,
819 M_PI/2., M_PI, "Total lapconf function seen by NS",
820 "Lapconf (theta=pi/2, phi=0)" ) ;
821
822 des_profile( confo_auto, 0., 10.,
823 M_PI/2., M_PI, "Self conformal factor of NS",
824 "Confo (theta=pi/2, phi=0)" ) ;
825
826 des_profile( confo_tot, 0., 10.,
827 M_PI/2., M_PI, "Total conformal factor seen by NS",
828 "Confo (theta=pi/2, phi=0)" ) ;
829
830 des_coupe_vect_z( shift_auto, 0., -3., 0.5, 3,
831 "Self shift vector of NS") ;
832
833 des_coupe_vect_z( shift_tot, 0., -3., 0.5, 3,
834 "Total shift vector seen by NS") ;
835 */
836 } // End of main loop
837
838 //=========================================================
839 // End of iteration
840 //=========================================================
841
842
843}
844}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition cmp.C:307
void annule_hard()
Sets the Cmp to zero in a hard way.
Definition cmp.C:341
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
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition scalar_pde.C:139
int get_dzpuis() const
Returns dzpuis.
Definition scalar.h:563
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition scalar.C:359
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Definition scalar.C:790
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
Definition scalar.C:397
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition scalar.h:560
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star ).
Definition star_bhns.h:192
Sym_tensor taij_auto
Part of the extrinsic curvature tensor generated by shift_auto , lapse_auto , and confo_auto .
Definition star_bhns.h:182
Scalar lapconf_auto
Lapconf function generated by the star.
Definition star_bhns.h:113
double velo_pot_bhns(const double &mass_bh, const double &sepa, bool kerrschild, int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
void equilibrium_bhns(double ent_c, const double &mass_bh, const double &sepa, bool kerrschild, int mer, int mermax_ns, int mermax_potvit, int mermax_poisson, int filter_r, int filter_r_s, int filter_p_s, double relax_poisson, double relax_potvit, double thres_adapt, double resize_ns, const Tbl &fact_resize, Tbl &diff)
Computes an equilibrium configuration.
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Definition star_bhns.h:220
Scalar confo_tot
Total conformal factor.
Definition star_bhns.h:163
Vector d_lapconf_auto
Derivative of the lapconf function generated by the star .
Definition star_bhns.h:130
Vector d_confo_auto
Derivative of the conformal factor generated by the star .
Definition star_bhns.h:168
double radius_p(double phi)
Radius of the star to the direction of and .
Scalar gam0
Lorentz factor between the co-orbiting observer and the Eulerian one.
Definition star_bhns.h:107
void hydro_euler_bhns(bool kerrschild, const double &mass_bh, const double &sepa)
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Scalar gam
Lorentz factor between the fluid and the co-orbiting observer.
Definition star_bhns.h:102
Scalar ssjm1_lapconf
Effective source at the previous step for the resolution of the Poisson equation for lapconf_auto .
Definition star_bhns.h:197
bool irrotational
true for an irrotational star, false for a corotating one
Definition star_bhns.h:72
Scalar pot_centri
Centrifugal potential.
Definition star_bhns.h:110
double chi_rp(double radius, double phi)
Sensitive indicator of the mass-shedding to the direction of , , .
Scalar lapconf_comp
Lapconf function generated by the companion black hole.
Definition star_bhns.h:116
Scalar taij_quad_auto
Part of the scalar generated by .
Definition star_bhns.h:187
double phi_min()
Azimuthal angle when the indicator of the mass-shedding takes its minimum chi_min.
Vector shift_auto
Shift vector generated by the star.
Definition star_bhns.h:138
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition star_bhns.h:210
Scalar confo_auto
Conformal factor generated by the star.
Definition star_bhns.h:157
Scalar ssjm1_confo
Effective source at the previous step for the resolution of the Poisson equation for confo_auto .
Definition star_bhns.h:202
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition star_bhns.h:93
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition star.h:198
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition star.C:465
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition star.h:201
Scalar press
Fluid pressure.
Definition star.h:194
double ray_eq_pi() const
Coordinate radius at , [r_unit].
Scalar ent
Log-enthalpy.
Definition star.h:190
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition star.h:207
Map & mp
Mapping associated with the star.
Definition star.h:180
int nzet
Number of domains of *mp occupied by the star.
Definition star.h:183
double ray_pole() const
Coordinate radius at [r_unit].
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 set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1176
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
Tensor field of valence 1.
Definition vector.h:188
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition vector.C:319
Scalar & set(int)
Read/write access to a component.
Definition vector.C:299
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Cmp exp(const Cmp &)
Exponential.
Definition cmp_math.C:273
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 pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
Cmp abs(const Cmp &)
Absolute value.
Definition cmp_math.C:413
Cmp log(const Cmp &)
Neperian logarithm.
Definition cmp_math.C:299
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition tensor.C:817
virtual void annule(int l_min, int l_max)
Sets the Tensor to zero in several domains.
Definition tensor.C:680
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition tensor.C:490
Lorene prototypes.
Definition app_hor.h:67
Standard units of space, time and mass.