LORENE
et_bin_bhns_extr_upmetr.C
1/*
2 * Methods Et_bin_bhns_extr::update_metric_extr_ks
3 * and Et_bin_bhns_extr::update_metric_extr_cf
4 *
5 * (see file et_bin_bhns_extr.h for documentation).
6 *
7 */
8
9/*
10 * Copyright (c) 2004-2005 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: et_bin_bhns_extr_upmetr.C,v 1.5 2016/12/05 16:17:52 j_novak Exp $
33 * $Log: et_bin_bhns_extr_upmetr.C,v $
34 * Revision 1.5 2016/12/05 16:17:52 j_novak
35 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36 *
37 * Revision 1.4 2014/10/13 08:52:55 j_novak
38 * Lorene classes and functions now belong to the namespace Lorene.
39 *
40 * Revision 1.3 2014/10/06 15:13:08 j_novak
41 * Modified #include directives to use c++ syntax.
42 *
43 * Revision 1.2 2005/02/28 23:16:37 k_taniguchi
44 * Modification to include the case of the conformally flat background metric
45 *
46 * Revision 1.1 2004/11/30 20:51:32 k_taniguchi
47 * *** empty log message ***
48 *
49 *
50 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_upmetr.C,v 1.5 2016/12/05 16:17:52 j_novak Exp $
51 *
52 */
53
54// C headers
55#include <cmath>
56
57// Lorene headers
58#include "et_bin_bhns_extr.h"
59#include "etoile.h"
60#include "coord.h"
61#include "unites.h"
62
63 //-----------------------------------------------------------//
64 // No relaxation for a fixed BH background //
65 //-----------------------------------------------------------//
66
67namespace Lorene {
69 const double& sepa)
70{
71
72 using namespace Unites ;
73
74 if (kerrschild) {
75
76 // Computation of quantities coming from the companion (K-S BH)
77 // ------------------------------------------------------------
78
79 const Coord& xx = mp.x ;
80 const Coord& yy = mp.y ;
81 const Coord& zz = mp.z ;
82
83 Tenseur r_bh(mp) ;
84 r_bh.set_etat_qcq() ;
85 r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
86 r_bh.set_std_base() ;
87
88 Tenseur xx_con(mp, 1, CON, ref_triad) ;
89 xx_con.set_etat_qcq() ;
90 xx_con.set(0) = xx + sepa ;
91 xx_con.set(1) = yy ;
92 xx_con.set(2) = zz ;
93 xx_con.set_std_base() ;
94
95 Tenseur xsr_con(mp, 1, CON, ref_triad) ;
96 xsr_con = xx_con / r_bh ;
97 xsr_con.set_std_base() ;
98
99 Tenseur msr(mp) ;
100 msr = ggrav * mass / r_bh ;
101 msr.set_std_base() ;
102
103 Tenseur lapse_bh(mp) ;
104 lapse_bh = 1. / sqrt( 1.+2.*msr ) ;
105 lapse_bh.set_std_base() ;
106
107 logn_comp.set_etat_qcq() ;
108 logn_comp.set() = log( lapse_bh() ) ;
109 logn_comp.set_std_base() ;
110
111 beta_comp.set_etat_qcq() ;
112 beta_comp.set() = log( lapse_bh() ) ;
113 // conformal factor of KS-BH is unity
114 beta_comp.set_std_base() ;
115
116 shift_comp.set_etat_qcq() ;
117
118 shift_comp.set(0) = -2.*lapse_bh()*lapse_bh()*msr()*xsr_con(0) ;
119 shift_comp.set(1) = -2.*lapse_bh()*lapse_bh()*msr()*xsr_con(1) ;
120 shift_comp.set(2) = -2.*lapse_bh()*lapse_bh()*msr()*xsr_con(2) ;
121
122 shift_comp.set_std_base() ;
123 shift_comp.set_triad( ref_triad ) ;
124
125 // Lapse function N
126 // ----------------
127
128 nnn = exp( unsurc2 * logn_auto ) * lapse_bh ;
129
130 nnn.set_std_base() ;
131
132 // Conformal factor A^2
133 // --------------------
134
135 a_car = exp ( 2.*unsurc2*(beta_auto - logn_auto) ) ;
136
137 a_car.set_std_base() ;
138
139 // Shift vector N^i
140 // ----------------
141
143
144 // Derivative of metric coefficients
145 // ----------------------------------
146
147 // ... (d/dX,d/dY,d/dZ)(logn_auto) :
148 d_logn_auto_regu = logn_auto_regu.gradient() ; // (d/dx, d/dy, d/dz)
149 d_logn_auto_regu.change_triad(ref_triad) ; // --> (d/dX, d/dY, d/dZ)
150
151 if ( *(d_logn_auto_div.get_triad()) != ref_triad ) {
152
153 // Change the basis from spherical coordinate to Cartesian one
154 d_logn_auto_div.change_triad( mp.get_bvect_cart() ) ;
155
156 // Change the basis from mapping coordinate to absolute one
157 d_logn_auto_div.change_triad( ref_triad ) ;
158
159 }
160
162
163 // ... (d/dX,d/dY,d/dZ)(beta_auto) :
164 d_beta_auto = beta_auto.gradient() ; // (d/dx, d/dy, d/dz)
165 d_beta_auto.change_triad(ref_triad) ; // --> (d/dX, d/dY, d/dZ)
166
167 if (relativistic) {
168 // ... extrinsic curvature (tkij_auto and akcar_auto)
169 extrinsic_curv_extr(mass, sepa) ;
170 }
171
172 // The derived quantities are obsolete
173 // -----------------------------------
174
176
177 }
178 else {
179
180 // Computation of quantities coming from the companion (CF Sch. BH)
181 // ----------------------------------------------------------------
182
183 const Coord& xx = mp.x ;
184 const Coord& yy = mp.y ;
185 const Coord& zz = mp.z ;
186
187 Tenseur r_bh(mp) ;
188 r_bh.set_etat_qcq() ;
189 r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
190 r_bh.set_std_base() ;
191
192 Tenseur msr(mp) ;
193 msr = ggrav * mass / r_bh ;
194 msr.set_std_base() ;
195
196 Tenseur lapse_bh(mp) ;
197 lapse_bh = (1.-0.5*msr) / (1.+0.5*msr) ;
198 lapse_bh.set_std_base() ;
199
200 logn_comp.set_etat_qcq() ;
201 logn_comp.set() = log( lapse_bh() ) ;
202 logn_comp.set_std_base() ;
203
204 Tenseur lappsi(mp) ;
205 lappsi = 1. - 0.25*msr*msr ;
206 lappsi.set_std_base() ;
207
208 beta_comp.set_etat_qcq() ;
209 beta_comp.set() = log( lappsi() ) ;
210 beta_comp.set_std_base() ;
211
212 shift_comp.set_etat_qcq() ;
213
214 shift_comp.set(0) = 0. ;
215 shift_comp.set(1) = 0. ;
216 shift_comp.set(2) = 0. ;
217
218 shift_comp.set_std_base() ;
219 shift_comp.set_triad( ref_triad ) ;
220
221 // Lapse function N
222 // ----------------
223
224 nnn = exp( unsurc2 * logn_auto ) * lapse_bh ;
225
226 nnn.set_std_base() ;
227
228 // Conformal factor A^2
229 // --------------------
230
232 - logn_auto - logn_comp) ) ;
233
234 a_car.set_std_base() ;
235
236 // Shift vector N^i
237 // ----------------
238
240
241 // Derivative of metric coefficients
242 // ----------------------------------
243
244 // ... (d/dX,d/dY,d/dZ)(logn_auto) :
245 d_logn_auto_regu = logn_auto_regu.gradient() ; // (d/dx, d/dy, d/dz)
246 d_logn_auto_regu.change_triad(ref_triad) ; // --> (d/dX, d/dY, d/dZ)
247
248 if ( *(d_logn_auto_div.get_triad()) != ref_triad ) {
249
250 // Change the basis from spherical coordinate to Cartesian one
251 d_logn_auto_div.change_triad( mp.get_bvect_cart() ) ;
252
253 // Change the basis from mapping coordinate to absolute one
254 d_logn_auto_div.change_triad( ref_triad ) ;
255
256 }
257
259
260 // ... (d/dX,d/dY,d/dZ)(beta_auto) :
261 d_beta_auto = beta_auto.gradient() ; // (d/dx, d/dy, d/dz)
262 d_beta_auto.change_triad(ref_triad) ; // --> (d/dX, d/dY, d/dZ)
263
264 if (relativistic) {
265 // ... extrinsic curvature (tkij_auto and akcar_auto)
266 extrinsic_curv_extr(mass, sepa) ;
267 }
268
269 // The derived quantities are obsolete
270 // -----------------------------------
271
273
274 }
275
276}
277}
Active physical coordinates and mapping derivatives.
Definition coord.h:90
bool kerrschild
Indicator of the background metric: true for the Kerr-Shild metric, false for the conformally flat on...
void update_metric_extr(const double &mass, const double &sepa)
Computes metric coefficients from known potentials, when the companion is a black hole with the Kerr-...
void extrinsic_curv_extr(const double &mass, const double &sepa)
Computes tkij_auto and akcar_auto from shift_auto , nnn and a_car .
Tenseur shift_comp
Part of the shift vector generated principaly by the companion star.
Definition etoile.h:898
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
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
Tenseur logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principaly by ...
Definition etoile.h:857
Tenseur beta_comp
Part of the logarithm of AN generated principaly by the companion star.
Definition etoile.h:877
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile_bin.C:450
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition etoile.h:892
Tenseur d_logn_auto_regu
Gradient of logn_auto_regu (Cartesian components with respect to ref_triad ).
Definition etoile.h:867
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
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
Map & mp
Mapping associated with the star.
Definition etoile.h:432
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 ).
Definition etoile.h:504
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:440
Tenseur shift
Total shift vector.
Definition etoile.h:515
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 unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition etoile.h:445
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
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Cmp exp(const Cmp &)
Exponential.
Definition cmp_math.C:273
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
Cmp log(const Cmp &)
Neperian logarithm.
Definition cmp_math.C:299
Lorene prototypes.
Definition app_hor.h:67
Standard units of space, time and mass.