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