LORENE
hole_bhns_upmetr.C
1 /*
2  * Method of class Hole_bhns to compute metric quantities from
3  * the companion neutron-star and total metric quantities
4  *
5  * (see file hole_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: hole_bhns_upmetr.C,v 1.4 2016/12/05 16:17:55 j_novak Exp $
33  * $Log: hole_bhns_upmetr.C,v $
34  * Revision 1.4 2016/12/05 16:17:55 j_novak
35  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36  *
37  * Revision 1.3 2014/10/13 08:53:00 j_novak
38  * Lorene classes and functions now belong to the namespace Lorene.
39  *
40  * Revision 1.2 2008/05/15 19:08:15 k_taniguchi
41  * Change of some parameters.
42  *
43  * Revision 1.1 2007/06/22 01:25:31 k_taniguchi
44  * *** empty log message ***
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_upmetr.C,v 1.4 2016/12/05 16:17:55 j_novak Exp $
48  *
49  */
50 
51 // C++ headers
52 //#include <>
53 
54 // C headers
55 //#include <>
56 
57 // Lorene headers
58 #include "hole_bhns.h"
59 #include "star_bhns.h"
60 #include "utilitaires.h"
61 #include "unites.h"
62 
63 namespace Lorene {
65  const Hole_bhns& hole_prev,
66  double relax) {
67 
68  // Fundamental constants and units
69  // -------------------------------
70  using namespace Unites ;
71 
72  //-----------------------------------------------------
73  // Computation of quantities coming from the companion
74  //-----------------------------------------------------
75 
76  const Map& mp_ns (star.get_mp()) ;
77 
78  double mass = ggrav * mass_bh ;
79 
80  // Lapconf function
81  // ----------------
82 
83  if ( (star.get_lapconf_auto()).get_etat() == ETATZERO ) {
85  }
86  else {
90  }
91 
92  // Shift vector
93  // ------------
94 
95  if ( (star.get_shift_auto())(2).get_etat() == ETATZERO ) {
96  assert( (star.get_shift_auto())(1).get_etat() == ETATZERO ) ;
97  assert( (star.get_shift_auto())(3).get_etat() == ETATZERO ) ;
98 
100  }
101  else {
104 
105  Vector comp_shift(star.get_shift_auto()) ;
106  comp_shift.change_triad(mp_ns.get_bvect_cart()) ;
107  comp_shift.change_triad(mp.get_bvect_cart()) ;
108 
109  assert( *(shift_comp.get_triad()) == *(comp_shift.get_triad()) ) ;
110 
111  (shift_comp.set(1)).import( comp_shift(1) ) ;
112  (shift_comp.set(2)).import( comp_shift(2) ) ;
113  (shift_comp.set(3)).import( comp_shift(3) ) ;
114 
116  }
117 
118  // Conformal factor
119  // ----------------
120 
121  if ( (star.get_confo_auto()).get_etat() == ETATZERO ) {
123  }
124  else {
126  confo_comp.import( star.get_confo_auto() ) ;
128  }
129 
130  //----------------------------------------------------
131  // Relaxation on lapconf_comp, shift_comp, confo_comp
132  //----------------------------------------------------
133 
134  double relax_jm1 = 1. - relax ;
135 
136  lapconf_comp = relax * lapconf_comp
137  + relax_jm1 * (hole_prev.lapconf_comp) ;
138 
139  shift_comp = relax * shift_comp + relax_jm1 * (hole_prev.shift_comp) ;
140 
141  confo_comp = relax * confo_comp + relax_jm1 * (hole_prev.confo_comp) ;
142 
143 
144  // Coordinates
145  // -----------
146  Scalar rr(mp) ;
147  rr = mp.r ;
148  rr.std_spectral_base() ;
149  Scalar st(mp) ;
150  st = mp.sint ;
151  st.std_spectral_base() ;
152  Scalar ct(mp) ;
153  ct = mp.cost ;
154  ct.std_spectral_base() ;
155  Scalar sp(mp) ;
156  sp = mp.sinp ;
157  sp.std_spectral_base() ;
158  Scalar cp(mp) ;
159  cp = mp.cosp ;
160  cp.std_spectral_base() ;
161 
162  Vector ll(mp, CON, mp.get_bvect_cart()) ;
163  ll.set_etat_qcq() ;
164  ll.set(1) = st % cp ;
165  ll.set(2) = st % sp ;
166  ll.set(3) = ct ;
167  ll.std_spectral_base() ;
168 
169 
170  if (kerrschild) {
171 
172  //----------------------------------------------
173  // Metric quantities from the analytic solution
174  //----------------------------------------------
175 
176  lapconf_auto_bh = 1. / sqrt(1.+2.*mass/rr) ;
180 
181  confo_auto_bh = 1. ;
183 
184  shift_auto_bh = 2. * lapconf_auto_bh*lapconf_auto_bh*mass * ll / rr ;
187 
188  //---------------------------------
189  // Derivative of metric quantities
190  //---------------------------------
191 
193  for (int i=1; i<=3; i++)
195 
197 
199  for (int i=1; i<=3; i++) {
200  d_lapconf_auto_bh.set(i) = pow(lapconf_auto_bh,3.) * mass * ll(i)
201  / rr / rr ;
202  }
206 
209 
211  for (int i=1; i<=3; i++) {
212  for (int j=1; j<=3; j++) {
213  d_shift_auto_rs.set(i,j) = shift_auto_rs(j).deriv(i) ;
214  }
215  }
216 
218 
220  for (int i=1; i<=3; i++) {
221  for (int j=1; j<=3; j++) {
223  *lapconf_auto_bh*mass
224  * (flat.con()(i,j)
225  - 2.*lapconf_auto_bh*lapconf_auto_bh*(1.+mass/rr)
226  * ll(i) * ll(j))
227  / rr / rr ;
228  }
229  }
233 
236 
238  for (int i=1; i<=3; i++)
240 
242 
244  for (int i=1; i<=3; i++)
245  d_confo_auto_bh.set(i) = 0. ;
246 
248 
251 
252  }
253  else { // Isotropic coordinates with the maximal slicing
254 
255  // Sets C/M^2 for each case of the lapse boundary condition
256  // --------------------------------------------------------
257  double cc ;
258 
259  if (bc_lapconf_nd) { // Neumann boundary condition
260  if (bc_lapconf_fs) { // First condition
261  // d(\alpha \psi)/dr = 0
262  // ---------------------
263  cc = 2. * (sqrt(13.) - 1.) / 3. ;
264  }
265  else { // Second condition
266  // d(\alpha \psi)/dr = (\alpha \psi)/(2 rah)
267  // -----------------------------------------
268  cc = 4. / 3. ;
269  }
270  }
271  else { // Dirichlet boundary condition
272  if (bc_lapconf_fs) { // First condition
273  // (\alpha \psi) = 1/2
274  // -------------------
275  cout << "!!!!! WARNING: Not yet prepared !!!!!" << endl ;
276  abort() ;
277  }
278  else { // Second condition
279  // (\alpha \psi) = 1/sqrt(2.) \psi_KS
280  // -----------------------------------
281  cout << "!!!!! WARNING: Not yet prepared !!!!!" << endl ;
282  abort() ;
283  // cc = 2. * sqrt(2.) ;
284  }
285  }
286 
287  //----------------------------------------------
288  // Metric quantities from the analytic solution
289  //----------------------------------------------
290 
291  Scalar r_are(mp) ;
293  r_are.std_spectral_base() ;
294 
295  lapconf_auto_bh = sqrt(1. - 2.*mass/r_are/rr
296  + cc*cc*pow(mass/r_are/rr, 4.)) * sqrt(r_are) ;
300 
301  confo_auto_bh = sqrt(r_are) ;
305 
306  shift_auto_bh = mass * mass * cc * ll / rr / rr / pow(r_are, 3.) ;
309  for (int i=1; i<=3; i++)
310  shift_auto_bh.set(i).raccord(1) ;
311 
312  //---------------------------------
313  // Derivative of metric quantities
314  //---------------------------------
315 
317  for (int i=1; i<=3; i++)
319 
321 
323  for (int i=1; i<=3; i++) {
324  d_lapconf_auto_bh.set(i) = sqrt(r_are)
325  * (mass/r_are/rr - 2.*cc*cc*pow(mass/r_are/rr,4.))
326  * ll(i) / rr
327  + 0.5 * sqrt(r_are)
328  * (sqrt(1. - 2.*mass/r_are/rr + cc*cc*pow(mass/r_are/rr,4.))-1.)
329  * sqrt(1. - 2.*mass/r_are/rr + cc*cc*pow(mass/r_are/rr,4.))
330  * ll(i) / rr ;
331  }
335 
339  for (int i=1; i<=3; i++)
340  d_lapconf_auto.set(i).raccord(1) ;
341 
343  for (int i=1; i<=3; i++) {
344  for (int j=1; j<=3; j++) {
345  d_shift_auto_rs.set(i,j) = shift_auto_rs(j).deriv(i) ;
346  }
347  }
348 
350 
352  for (int i=1; i<=3; i++) {
353  for (int j=1; j<=3; j++) {
354  d_shift_auto_bh.set(i,j) =
355  mass*mass*cc*(flat.con()(i,j)
356  -3.*sqrt(1. - 2.*mass/r_are/rr
357  +cc*cc*pow(mass/r_are/rr,4.))
358  *ll(i)*ll(j))
359  / pow(r_are*rr,3.) ;
360  }
361  }
365 
369  for (int i=1; i<=3; i++) {
370  for (int j=1; j<=3; j++) {
371  d_shift_auto.set(i,j).raccord(1) ;
372  }
373  }
374 
376  for (int i=1; i<=3; i++)
378 
380 
382  for (int i=1; i<=3; i++) {
383  d_confo_auto_bh.set(i) = 0.5*sqrt(r_are)
384  *(sqrt(1. - 2.*mass/r_are/rr +cc*cc*pow(mass/r_are/rr,4.)) - 1.)
385  * ll(i) / rr ;
386  }
390 
394  for (int i=1; i<=3; i++)
395  d_confo_auto.set(i).raccord(1) ;
396 
397  }
398 
399  //-------------------------
400  // Total metric quantities
401  //-------------------------
402 
406  lapconf_auto.raccord(1) ;
407 
411  lapconf_tot.raccord(1) ;
412 
416  for (int i=1; i<=3; i++) {
417  shift_auto.set(i).raccord(1) ;
418  }
419 
423  for (int i=1; i<=3; i++) {
424  shift_tot.set(i).raccord(1) ;
425  }
426 
430  confo_auto.raccord(1) ;
431 
435  confo_tot.raccord(1) ;
436 
439 
442 
443  // The derived quantities are obsolete
444  // -----------------------------------
445 
446  del_deriv() ;
447 
448 }
449 }
Map & mp
Mapping associated with the black hole.
Definition: blackhole.h:80
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:675
Scalar lapconf_auto_bh
Part of the lapconf function from the analytic background.
Definition: hole_bhns.h:92
Class for stars in black hole-neutron star binaries.
Definition: star_bhns.h:59
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition: tensor.C:490
Vector d_confo_auto
Derivative of the conformal factor generated by the black hole.
Definition: hole_bhns.h:182
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
Scalar confo_auto
Conformal factor generated by the black hole.
Definition: hole_bhns.h:163
Lorene prototypes.
Definition: app_hor.h:67
double mass_bh
Gravitational mass of BH.
Definition: blackhole.h:88
Standard units of space, time and mass.
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
Vector d_lapconf_auto_bh
Derivative of the part of the lapconf function from the analytic background.
Definition: hole_bhns.h:117
bool kerrschild
true for a Kerr-Schild background, false for a conformally flat background
Definition: blackhole.h:85
Base class for coordinate mappings.
Definition: map.h:688
Vector shift_auto
Shift vector generated by the black hole.
Definition: hole_bhns.h:132
Vector d_lapconf_auto
Derivative of the lapconf function generated by the black hole.
Definition: hole_bhns.h:120
Tensor d_shift_auto
Derivative of the shift vector generated by the black hole.
Definition: hole_bhns.h:151
Vector shift_tot
Total shift vector ;.
Definition: hole_bhns.h:138
const Scalar & get_lapconf_auto() const
Returns the part of the lapconf function generated by the star.
Definition: star_bhns.h:339
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
const Map & get_mp() const
Returns the mapping.
Definition: star.h:355
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
Tensor field of valence 1.
Definition: vector.h:188
Scalar lapconf_comp
Lapconf function generated by the companion star.
Definition: hole_bhns.h:98
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: vector.C:322
void raccord(int n)
Performs the matching of the nucleus with respect to the first shell.
Tensor d_shift_auto_bh
Derivative of the part of the shift vector from the analytic background.
Definition: hole_bhns.h:148
Coord sint
Definition: map.h:739
Scalar confo_auto_rs
Part of the conformal factor from the numerical computation.
Definition: hole_bhns.h:157
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition: metric_flat.C:156
const Vector & get_shift_auto() const
Returns the part of the shift vector generated by the star.
Definition: star_bhns.h:364
Scalar lapse_tot
Total lapse function.
Definition: hole_bhns.h:107
Scalar lapse_auto
Lapse function of the "black hole" part.
Definition: hole_bhns.h:104
Vector shift_auto_bh
Part of the shift vector from the analytic background.
Definition: hole_bhns.h:129
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tensor.h:879
Scalar confo_comp
Conformal factor generated by the companion star.
Definition: hole_bhns.h:166
Vector shift_auto_rs
Part of the shift vector from the numerical computation.
Definition: hole_bhns.h:126
Coord sinp
Definition: map.h:741
Vector shift_comp
Shift vector generated by the companion star.
Definition: hole_bhns.h:135
const Scalar & deriv(int i) const
Returns of *this , where .
Definition: scalar_deriv.C:359
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
const Scalar & get_confo_auto() const
Returns the part of the conformal factor generated by the star.
Definition: star_bhns.h:383
Vector d_lapconf_auto_rs
Derivative of the part of the lapconf function from the numerical computation.
Definition: hole_bhns.h:112
Tensor d_shift_auto_rs
Derivative of the part of the shift vector from the numerical computation.
Definition: hole_bhns.h:143
Vector d_confo_auto_bh
Derivative of the part of the conformal factor from the analytic background.
Definition: hole_bhns.h:179
Class for black holes in black hole-neutron star binaries.
Definition: hole_bhns.h:65
Scalar lapconf_auto
Lapconf function generated by the black hole.
Definition: hole_bhns.h:95
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:825
Scalar confo_tot
Total conformal factor.
Definition: hole_bhns.h:169
void import(const Scalar &ci)
Assignment to another Scalar defined on a different mapping.
Definition: scalar_import.C:71
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:809
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:663
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the black hol...
Definition: blackhole.h:143
Coord cosp
Definition: map.h:742
Scalar lapconf_auto_rs
Part of the lapconf function from the numerical computation.
Definition: hole_bhns.h:89
bool bc_lapconf_fs
true for the first type BC for the lapconf function, false for the second type BH ...
Definition: hole_bhns.h:78
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tensor.C:528
Vector d_confo_auto_rs
Derivative of the part of the conformal factor from the numerical computation.
Definition: hole_bhns.h:174
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:506
bool bc_lapconf_nd
true for the Neumann type BC for the lapconf function, false for the Dirichlet type BH ...
Definition: hole_bhns.h:73
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
const Scalar r_coord(bool neumann, bool first) const
Expresses the areal radial coordinate by that in spatially isotropic coordinates. ...
void update_metric_bhns(const Star_bhns &star, const Hole_bhns &hole_prev, double relax)
Computes metric coefficients from known potentials with relaxation when the companion is a black hole...
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: tensor.C:935
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: hole_bhns.C:395
Scalar lapconf_tot
Total lapconf function.
Definition: hole_bhns.h:101
Coord r
r coordinate centered on the grid
Definition: map.h:736
Scalar confo_auto_bh
Part of the conformal factor from the analytic background.
Definition: hole_bhns.h:160
Coord cost
Definition: map.h:740