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 char et_bin_bhns_extr_um_der_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_um_der.C,v 1.4 2014/10/13 08:52:55 j_novak Exp $" ;
30 
31 /*
32  * $Id: et_bin_bhns_extr_um_der.C,v 1.4 2014/10/13 08:52:55 j_novak Exp $
33  * $Log: et_bin_bhns_extr_um_der.C,v $
34  * Revision 1.4 2014/10/13 08:52:55 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.3 2014/10/06 15:13:08 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.2 2005/02/28 23:15:52 k_taniguchi
41  * Modification to include the case of the conformally flat background metric
42  *
43  * Revision 1.1 2004/11/30 20:51:09 k_taniguchi
44  * *** empty log message ***
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_um_der.C,v 1.4 2014/10/13 08:52:55 j_novak Exp $
48  *
49  */
50 
51 // C headers
52 #include <cmath>
53 
54 // Lorene headers
55 #include "et_bin_bhns_extr.h"
56 #include "etoile.h"
57 #include "coord.h"
58 #include "unites.h"
59 
60 namespace Lorene {
62  const double& sepa) {
63 
64  using namespace Unites ;
65 
66  if (kerrschild) {
67 
68  // Computation of quantities coming from the companion (K-S BH)
69  // ------------------------------------------------------------
70 
71  const Coord& xx = mp.x ;
72  const Coord& yy = mp.y ;
73  const Coord& zz = mp.z ;
74 
75  Tenseur r_bh(mp) ;
76  r_bh.set_etat_qcq() ;
77  r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
78  r_bh.set_std_base() ;
79 
80  Tenseur xx_cov(mp, 1, COV, ref_triad) ;
81  xx_cov.set_etat_qcq() ;
82  xx_cov.set(0) = xx + sepa ;
83  xx_cov.set(1) = yy ;
84  xx_cov.set(2) = zz ;
85  xx_cov.set_std_base() ;
86 
87  Tenseur xsr_cov(mp, 1, COV, ref_triad) ;
88  xsr_cov = xx_cov / r_bh ;
89  xsr_cov.set_std_base() ;
90 
91  Tenseur xx_con(mp, 1, CON, ref_triad) ;
92  xx_con.set_etat_qcq() ;
93  xx_con.set(0) = xx + sepa ;
94  xx_con.set(1) = yy ;
95  xx_con.set(2) = zz ;
96  xx_con.set_std_base() ;
97 
98  Tenseur xsr_con(mp, 1, CON, ref_triad) ;
99  xsr_con = xx_con / r_bh ;
100  xsr_con.set_std_base() ;
101 
102  Tenseur msr(mp) ;
103  msr = ggrav * mass / r_bh ;
104  msr.set_std_base() ;
105 
106  // Computation of d_logn_comp
107  // --------------------------
108 
109  Tenseur lapse_bh2(mp) ; // lapse_bh * lapse_bh
110  lapse_bh2 = 1. / (1.+2.*msr) ;
111  lapse_bh2.set_std_base() ;
112 
114 
115  d_logn_comp.set(0) = lapse_bh2()%msr()%xsr_cov(0)/r_bh() ;
116  d_logn_comp.set(1) = lapse_bh2()%msr()%xsr_cov(1)/r_bh() ;
117  d_logn_comp.set(2) = lapse_bh2()%msr()%xsr_cov(2)/r_bh() ;
118 
121 
122  // Computation of d_beta_comp (Just the same that d_logn_comp)
123  // --------------------------
124 
127 
128  // Computation of tkij_comp
129  // ------------------------
134  // Components of shift_comp with respect to the Cartesian triad
135  // (d/dx, d/dy, d/dz) of the mapping :
136 
137  // Computation of A^2 A^{ij}
138  // -------------------------
140 
141  for (int i=0; i<3; i++) {
142  for (int j=i; j<3; j++) {
143  tkij_comp.set(i, j) = - 3.*xsr_con(i) % xsr_con(j)
144  - 2*msr() % xsr_con(i) % xsr_con(j) ;
145  }
146  tkij_comp.set(i, i) += 1. + 2.*msr() ;
147  }
148 
149  tkij_comp = (double(2)/double(3)) * pow(lapse_bh2, 3.) % msr
150  % (2.*tkij_comp +3.*msr % tkij_comp) / nnn / r_bh ;
151 
154 
155  if (relativistic) {
156  // Computation of akcar_comp
157  // -------------------------
158 
159  Tenseur lapse_bh8(mp) ;
160  lapse_bh8 = 1. / pow(1.+2.*msr, 4.) ;
161  lapse_bh8.set_std_base() ;
162 
164  akcar_comp.set() = 0 ;
165 
166  akcar_comp.set() = 8.*lapse_bh8()
167  * pow(2.*msr()+3.*msr()*msr(), 2.) / 3.
168  / nnn() / nnn() / r_bh() / r_bh() ;
169 
172 
173  }
174 
175  // The derived quantities are obsolete
176  // -----------------------------------
177 
179 
180  }
181  else {
182 
183  // Computation of quantities coming from the companion (CF Sch. BH)
184  // ----------------------------------------------------------------
185 
186  const Coord& xx = mp.x ;
187  const Coord& yy = mp.y ;
188  const Coord& zz = mp.z ;
189 
190  Tenseur r_bh(mp) ;
191  r_bh.set_etat_qcq() ;
192  r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
193  r_bh.set_std_base() ;
194 
195  Tenseur xx_cov(mp, 1, COV, ref_triad) ;
196  xx_cov.set_etat_qcq() ;
197  xx_cov.set(0) = xx + sepa ;
198  xx_cov.set(1) = yy ;
199  xx_cov.set(2) = zz ;
200  xx_cov.set_std_base() ;
201 
202  Tenseur msr(mp) ;
203  msr = ggrav * mass / r_bh ;
204  msr.set_std_base() ;
205 
206  // Computation of d_logn_comp
207  // --------------------------
208 
209  Tenseur tmp(mp) ;
210  tmp = 1. / (1. - 0.25*msr*msr) ;
211  tmp.set_std_base() ;
212 
214 
215  d_logn_comp.set(0) = tmp()%msr()%xx_cov(0)/r_bh()/r_bh() ;
216  d_logn_comp.set(1) = tmp()%msr()%xx_cov(1)/r_bh()/r_bh() ;
217  d_logn_comp.set(2) = tmp()%msr()%xx_cov(2)/r_bh()/r_bh() ;
218 
221 
222  // Computation of d_beta_comp
223  // --------------------------
224 
226 
227  d_beta_comp.set(0) = 0.5*tmp()%msr()%msr()%xx_cov(0)/r_bh()/r_bh() ;
228  d_beta_comp.set(1) = 0.5*tmp()%msr()%msr()%xx_cov(1)/r_bh()/r_bh() ;
229  d_beta_comp.set(2) = 0.5*tmp()%msr()%msr()%xx_cov(2)/r_bh()/r_bh() ;
230 
233 
234  // Computation of tkij_comp
235  // ------------------------
236 
237  // Computation of A^2 A^{ij}
238  // -------------------------
240 
241  for (int i=0; i<3; i++) {
242  for (int j=i; j<3; j++) {
243  tkij_comp.set(i, j) = 0. ;
244  }
245  }
246 
249 
250  if (relativistic) {
251  // Computation of akcar_comp
252  // -------------------------
253 
255  akcar_comp.set() = 0 ;
257 
258  }
259 
260  // The derived quantities are obsolete
261  // -----------------------------------
262 
264 
265  }
266 
267 }
268 }
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:828
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tenseur.C:674
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1170
Lorene prototypes.
Definition: app_hor.h:64
Standard units of space, time and mass.
Tenseur nnn
Total lapse function.
Definition: etoile.h:509
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:824
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:884
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
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_bin.C:447
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:944
Map & mp
Mapping associated with the star.
Definition: etoile.h:429
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:348
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
Tenseur a_car
Total conformal factor .
Definition: etoile.h:515
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:437
Coord y
y coordinate centered on the grid
Definition: map.h:727
Coord x
x coordinate centered on the grid
Definition: map.h:726
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:636
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:869
Coord z
z coordinate centered on the grid
Definition: map.h:728
Tenseur_sym tkij_comp
Part of the extrinsic curvature tensor generated by shift_comp .
Definition: etoile.h:932
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:301