LORENE
et_bin_bhns_extr_kinema.C
1 /*
2  * Method Et_bin_bhns_extr::kinematics_extr_ks
3  * and Et_bin_bhns_extr::kinematics_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_kinema.C,v 1.4 2016/12/05 16:17:52 j_novak Exp $
33  * $Log: et_bin_bhns_extr_kinema.C,v $
34  * Revision 1.4 2016/12/05 16:17:52 j_novak
35  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36  *
37  * Revision 1.3 2014/10/13 08:52:55 j_novak
38  * Lorene classes and functions now belong to the namespace Lorene.
39  *
40  * Revision 1.2 2005/02/28 23:15:09 k_taniguchi
41  * Modification to include the case of the conformally flat background metric
42  *
43  * Revision 1.1 2004/11/30 20:49:58 k_taniguchi
44  * *** empty log message ***
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_kinema.C,v 1.4 2016/12/05 16:17:52 j_novak Exp $
48  *
49  */
50 
51 // Lorene headers
52 #include "et_bin_bhns_extr.h"
53 #include "etoile.h"
54 #include "coord.h"
55 #include "unites.h"
56 
57 namespace Lorene {
58 void Et_bin_bhns_extr::kinematics_extr(double omega, const double& mass,
59  const double& sepa) {
60 
61  using namespace Unites ;
62 
63  if (kerrschild) {
64 
65  int nz = mp.get_mg()->get_nzone() ;
66  int nzm1 = nz - 1 ;
67 
68  // --------------------
69  // Computation of B^i/N
70  // --------------------
71 
72  // 1/ Computation of - omega m^i
73 
74  const Coord& xa = mp.xa ;
75  const Coord& ya = mp.ya ;
76 
77  bsn.set_etat_qcq() ;
78 
79  bsn.set(0) = omega * ya ;
80  bsn.set(1) = - omega * xa ;
81  bsn.set(2) = 0 ;
82 
83  bsn.annule(nzm1, nzm1) ; // set to zero in the ZEC
84 
85  // 2/ Addition of shift and division by lapse
86 
87  bsn = ( bsn + shift ) / nnn ;
88 
89  bsn.annule(nzm1, nzm1) ; // set to zero in the ZEC
90  bsn.set_std_base() ; // set the bases for spectral expansions
91 
92  //-------------------------
93  // Centrifugal potentatial
94  //-------------------------
95 
96  const Coord& xx = mp.x ;
97  const Coord& yy = mp.y ;
98  const Coord& zz = mp.z ;
99 
100  Tenseur r_bh(mp) ;
101  r_bh.set_etat_qcq() ;
102  r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
103  r_bh.set_std_base() ;
104 
105  Tenseur xx_cov(mp, 1, COV, ref_triad) ;
106  xx_cov.set_etat_qcq() ;
107  xx_cov.set(0) = xx + sepa ;
108  xx_cov.set(1) = yy ;
109  xx_cov.set(2) = zz ;
110  xx_cov.set_std_base() ;
111 
112  Tenseur xsr_cov(mp, 1, COV, ref_triad) ;
113  xsr_cov = xx_cov / r_bh ;
114  xsr_cov.set_std_base() ;
115 
116  Tenseur msr(mp) ;
117  msr = ggrav * mass / r_bh ;
118  msr.set_std_base() ;
119 
120  if (relativistic) {
121 
122  // Lorentz factor between the co-orbiting observer and
123  // the Eulerian one
124 
125  Tenseur tmp1(mp) ;
126  tmp1.set_etat_qcq() ;
127  tmp1.set() = 0 ;
128  tmp1.set_std_base() ;
129 
130  for (int i=0; i<3; i++)
131  tmp1.set() += xsr_cov(i) % bsn(i) ;
132 
133  tmp1.set_std_base() ;
134 
135  Tenseur tmp2 = 2.*msr % tmp1 % tmp1 ;
136  tmp2.set_std_base() ;
137 
138  for (int i=0; i<3; i++)
139  tmp2.set() += bsn(i) % bsn(i) ;
140 
141  tmp2 = a_car % tmp2 ;
142 
143  Tenseur gam0 = 1 / sqrt( 1 - tmp2 ) ;
144 
145  pot_centri = - log( gam0 ) ;
146 
147  }
148  else {
149  cout << "BH-NS binary system should be relativistic !!!" << endl ;
150  abort() ;
151  }
152 
153  pot_centri.annule(nzm1, nzm1) ; // set to zero in the external domain
154  pot_centri.set_std_base() ; // set the bases for spectral expansions
155 
156  // The derived quantities are obsolete
157  // -----------------------------------
158 
160 
161  }
162  else {
163 
164  int nz = mp.get_mg()->get_nzone() ;
165  int nzm1 = nz - 1 ;
166 
167  // --------------------
168  // Computation of B^i/N
169  // --------------------
170 
171  // 1/ Computation of - omega m^i
172 
173  const Coord& xa = mp.xa ;
174  const Coord& ya = mp.ya ;
175 
176  bsn.set_etat_qcq() ;
177 
178  bsn.set(0) = omega * ya ;
179  bsn.set(1) = - omega * xa ;
180  bsn.set(2) = 0 ;
181 
182  bsn.annule(nzm1, nzm1) ; // set to zero in the ZEC
183 
184  // 2/ Addition of shift and division by lapse
185 
186  bsn = ( bsn + shift ) / nnn ;
187 
188  bsn.annule(nzm1, nzm1) ; // set to zero in the ZEC
189  bsn.set_std_base() ; // set the bases for spectral expansions
190 
191  //-------------------------
192  // Centrifugal potentatial
193  //-------------------------
194 
195  if (relativistic) {
196 
197  // Lorentz factor between the co-orbiting observer and
198  // the Eulerian one
199 
200  Tenseur tmp(mp) ;
201  tmp.set_etat_qcq() ;
202  tmp.set() = 0. ;
203  tmp.set_std_base() ;
204 
205  for (int i=0; i<3; i++)
206  tmp.set() += bsn(i) % bsn(i) ;
207 
208  tmp = a_car % tmp ;
209 
210  Tenseur gam0 = 1 / sqrt( 1 - tmp ) ;
211 
212  pot_centri = - log( gam0 ) ;
213 
214  }
215  else {
216  cout << "BH-NS binary system should be relativistic !!!" << endl ;
217  abort() ;
218  }
219 
220  pot_centri.annule(nzm1, nzm1) ; // set to zero in the external domain
221  pot_centri.set_std_base() ; // set the bases for spectral expansions
222 
223  // The derived quantities are obsolete
224  // -----------------------------------
225 
227 
228  }
229 
230 }
231 }
Coord xa
Absolute x coordinate.
Definition: map.h:748
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
void annule(int l)
Sets the Tenseur to zero in a given domain.
Definition: tenseur.C:916
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
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
Tenseur pot_centri
Centrifugal potential.
Definition: etoile.h:956
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:783
Tenseur nnn
Total lapse function.
Definition: etoile.h:512
void kinematics_extr(double omega, const double &mass, const double &sepa)
Computes the quantities bsn and pot_centri in the Kerr-Schild background metric or in the conformally...
bool kerrschild
Indicator of the background metric: true for the Kerr-Shild metric, false for the conformally flat on...
Tenseur shift
Total shift vector.
Definition: etoile.h:515
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_bin.C:450
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition: etoile.h:953
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
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
Coord ya
Absolute y coordinate.
Definition: map.h:749
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
Coord z
z coordinate centered on the grid
Definition: map.h:746
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304