LORENE
et_bin_bhns_extr_um_der.C
1 /*
2  * Methods Et_bin_bhns_extr::update_metric_der_comp_extr_ks
3  * and Et_bin_bhns_extr::update_metric_der_comp_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_um_der.C,v 1.5 2016/12/05 16:17:52 j_novak Exp $
33  * $Log: et_bin_bhns_extr_um_der.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:15:52 k_taniguchi
44  * Modification to include the case of the conformally flat background metric
45  *
46  * Revision 1.1 2004/11/30 20:51:09 k_taniguchi
47  * *** empty log message ***
48  *
49  *
50  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_um_der.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 namespace Lorene {
65  const double& sepa) {
66 
67  using namespace Unites ;
68 
69  if (kerrschild) {
70 
71  // Computation of quantities coming from the companion (K-S BH)
72  // ------------------------------------------------------------
73 
74  const Coord& xx = mp.x ;
75  const Coord& yy = mp.y ;
76  const Coord& zz = mp.z ;
77 
78  Tenseur r_bh(mp) ;
79  r_bh.set_etat_qcq() ;
80  r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
81  r_bh.set_std_base() ;
82 
83  Tenseur xx_cov(mp, 1, COV, ref_triad) ;
84  xx_cov.set_etat_qcq() ;
85  xx_cov.set(0) = xx + sepa ;
86  xx_cov.set(1) = yy ;
87  xx_cov.set(2) = zz ;
88  xx_cov.set_std_base() ;
89 
90  Tenseur xsr_cov(mp, 1, COV, ref_triad) ;
91  xsr_cov = xx_cov / r_bh ;
92  xsr_cov.set_std_base() ;
93 
94  Tenseur xx_con(mp, 1, CON, ref_triad) ;
95  xx_con.set_etat_qcq() ;
96  xx_con.set(0) = xx + sepa ;
97  xx_con.set(1) = yy ;
98  xx_con.set(2) = zz ;
99  xx_con.set_std_base() ;
100 
101  Tenseur xsr_con(mp, 1, CON, ref_triad) ;
102  xsr_con = xx_con / r_bh ;
103  xsr_con.set_std_base() ;
104 
105  Tenseur msr(mp) ;
106  msr = ggrav * mass / r_bh ;
107  msr.set_std_base() ;
108 
109  // Computation of d_logn_comp
110  // --------------------------
111 
112  Tenseur lapse_bh2(mp) ; // lapse_bh * lapse_bh
113  lapse_bh2 = 1. / (1.+2.*msr) ;
114  lapse_bh2.set_std_base() ;
115 
117 
118  d_logn_comp.set(0) = lapse_bh2()%msr()%xsr_cov(0)/r_bh() ;
119  d_logn_comp.set(1) = lapse_bh2()%msr()%xsr_cov(1)/r_bh() ;
120  d_logn_comp.set(2) = lapse_bh2()%msr()%xsr_cov(2)/r_bh() ;
121 
124 
125  // Computation of d_beta_comp (Just the same that d_logn_comp)
126  // --------------------------
127 
130 
131  // Computation of tkij_comp
132  // ------------------------
137  // Components of shift_comp with respect to the Cartesian triad
138  // (d/dx, d/dy, d/dz) of the mapping :
139 
140  // Computation of A^2 A^{ij}
141  // -------------------------
143 
144  for (int i=0; i<3; i++) {
145  for (int j=i; j<3; j++) {
146  tkij_comp.set(i, j) = - 3.*xsr_con(i) % xsr_con(j)
147  - 2*msr() % xsr_con(i) % xsr_con(j) ;
148  }
149  tkij_comp.set(i, i) += 1. + 2.*msr() ;
150  }
151 
152  tkij_comp = (double(2)/double(3)) * pow(lapse_bh2, 3.) % msr
153  % (2.*tkij_comp +3.*msr % tkij_comp) / nnn / r_bh ;
154 
157 
158  if (relativistic) {
159  // Computation of akcar_comp
160  // -------------------------
161 
162  Tenseur lapse_bh8(mp) ;
163  lapse_bh8 = 1. / pow(1.+2.*msr, 4.) ;
164  lapse_bh8.set_std_base() ;
165 
167  akcar_comp.set() = 0 ;
168 
169  akcar_comp.set() = 8.*lapse_bh8()
170  * pow(2.*msr()+3.*msr()*msr(), 2.) / 3.
171  / nnn() / nnn() / r_bh() / r_bh() ;
172 
175 
176  }
177 
178  // The derived quantities are obsolete
179  // -----------------------------------
180 
182 
183  }
184  else {
185 
186  // Computation of quantities coming from the companion (CF Sch. BH)
187  // ----------------------------------------------------------------
188 
189  const Coord& xx = mp.x ;
190  const Coord& yy = mp.y ;
191  const Coord& zz = mp.z ;
192 
193  Tenseur r_bh(mp) ;
194  r_bh.set_etat_qcq() ;
195  r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
196  r_bh.set_std_base() ;
197 
198  Tenseur xx_cov(mp, 1, COV, ref_triad) ;
199  xx_cov.set_etat_qcq() ;
200  xx_cov.set(0) = xx + sepa ;
201  xx_cov.set(1) = yy ;
202  xx_cov.set(2) = zz ;
203  xx_cov.set_std_base() ;
204 
205  Tenseur msr(mp) ;
206  msr = ggrav * mass / r_bh ;
207  msr.set_std_base() ;
208 
209  // Computation of d_logn_comp
210  // --------------------------
211 
212  Tenseur tmp(mp) ;
213  tmp = 1. / (1. - 0.25*msr*msr) ;
214  tmp.set_std_base() ;
215 
217 
218  d_logn_comp.set(0) = tmp()%msr()%xx_cov(0)/r_bh()/r_bh() ;
219  d_logn_comp.set(1) = tmp()%msr()%xx_cov(1)/r_bh()/r_bh() ;
220  d_logn_comp.set(2) = tmp()%msr()%xx_cov(2)/r_bh()/r_bh() ;
221 
224 
225  // Computation of d_beta_comp
226  // --------------------------
227 
229 
230  d_beta_comp.set(0) = 0.5*tmp()%msr()%msr()%xx_cov(0)/r_bh()/r_bh() ;
231  d_beta_comp.set(1) = 0.5*tmp()%msr()%msr()%xx_cov(1)/r_bh()/r_bh() ;
232  d_beta_comp.set(2) = 0.5*tmp()%msr()%msr()%xx_cov(2)/r_bh()/r_bh() ;
233 
236 
237  // Computation of tkij_comp
238  // ------------------------
239 
240  // Computation of A^2 A^{ij}
241  // -------------------------
243 
244  for (int i=0; i<3; i++) {
245  for (int j=i; j<3; j++) {
246  tkij_comp.set(i, j) = 0. ;
247  }
248  }
249 
252 
253  if (relativistic) {
254  // Computation of akcar_comp
255  // -------------------------
256 
258  akcar_comp.set() = 0 ;
260 
261  }
262 
263  // The derived quantities are obsolete
264  // -----------------------------------
265 
267 
268  }
269 
270 }
271 }
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur &#39;s...
Definition: etoile.h:831
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tenseur.C:690
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
Tenseur nnn
Total lapse function.
Definition: etoile.h:512
bool kerrschild
Indicator of the background metric: true for the Kerr-Shild metric, false for the conformally flat on...
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:887
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:684
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_bin.C:450
void update_metric_der_comp_extr(const double &mass, const double &sepa)
Computes the derivative of metric functions related to the companion black hole with the Kerr-Schild ...
Tenseur akcar_comp
Part of the scalar generated by shift_auto and shift_comp , i.e.
Definition: etoile.h:947
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
Tenseur a_car
Total conformal factor .
Definition: etoile.h:518
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:440
Coord y
y coordinate centered on the grid
Definition: map.h:745
Coord x
x coordinate centered on the grid
Definition: map.h:744
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:872
Coord z
z coordinate centered on the grid
Definition: map.h:746
Tenseur_sym tkij_comp
Part of the extrinsic curvature tensor generated by shift_comp .
Definition: etoile.h:935
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304