LORENE
et_bin_nsbh_upmetr.C
1/*
2 * Methods Et_bin_nsbh::update_metric
3 *
4 * (see file et_bin_nsbh.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2003 Philippe Grandclement
10 * 2003 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_nsbh_upmetr.C,v 1.5 2016/12/05 16:17:53 j_novak Exp $
33 * $Log: et_bin_nsbh_upmetr.C,v $
34 * Revision 1.5 2016/12/05 16:17:53 j_novak
35 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36 *
37 * Revision 1.4 2014/10/13 08:52:56 j_novak
38 * Lorene classes and functions now belong to the namespace Lorene.
39 *
40 * Revision 1.3 2007/04/24 20:14:45 f_limousin
41 * Implementation of Dirichlet and Neumann BC for the lapse
42 *
43 * Revision 1.2 2005/08/29 15:10:17 p_grandclement
44 * Addition of things needed :
45 * 1) For BBH with different masses
46 * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
47 * WORKING YET !!!
48 *
49 * Revision 1.1 2003/10/24 12:28:39 k_taniguchi
50 * Method of update metric for the BH companion
51 *
52 *
53 *
54 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_upmetr.C,v 1.5 2016/12/05 16:17:53 j_novak Exp $
55 *
56 */
57
58// Lorene headers
59#include "et_bin_nsbh.h"
60#include "bhole.h"
61
62 //----------------------------------//
63 // Version a BH companion //
64 //----------------------------------//
65
66namespace Lorene {
68
69 // Computation of quantities coming from the companion
70 // ---------------------------------------------------
71
72 // Computes N_comp ---> stored in logn_comp
73 if ( (comp.get_n_auto()).get_etat() == ETATZERO ) {
74 n_comp.set_etat_zero() ;
75 }
76 else{
77 n_comp.set_etat_qcq() ;
78 (n_comp.set()).import( comp.get_n_auto()() ) ;
79 n_comp.set_std_base() ; // set the bases for spectral expansions
80 }
81
82
83 // Computes Psi_comp ---> stored in beta_comp
84 if ( (comp.get_psi_auto()).get_etat() == ETATZERO ) {
85 confpsi_comp.set_etat_zero() ;
86 }
87 else{
88 confpsi_comp.set_etat_qcq() ;
89 (confpsi_comp.set()).import( comp.get_psi_auto()() ) ;
90 confpsi_comp.set_std_base() ; // set the bases for spectral expansions
91 }
92
93
94 // Computes N^i_comp ---> stored in shift_comp
95 if ( (comp.get_shift_auto()).get_etat() == ETATZERO ) {
96 shift_comp.set_etat_zero() ;
97 }
98 else{
99 shift_comp.set_etat_qcq() ;
100
101 (shift_comp.set(0)).import( comp.get_shift_auto()(0) ) ; // N^x antisym
102 (shift_comp.set(1)).import( comp.get_shift_auto()(1) ) ; // N^y sym.
103 (shift_comp.set(2)).import( comp.get_shift_auto()(2) ) ; // N^z anisym
104
105 shift_comp.set_std_base() ; // set the bases for spectral expansions
106 }
107
108 // Lapse function N
109 // ----------------
110
111 nnn = n_auto + n_comp ;
112
113 // Conformal factor confpsi
114 // ------------------------
115
117
118 confpsi.set_std_base() ; // set the bases for spectral expansions
119
120 // Conformal factor A^2
121 // ---------------------
122
123 a_car = pow(confpsi, 4.);
124 a_car.set_std_base() ; // set the bases for spectral expansions
125
126 // Shift vector N^i
127 // ----------------
128
130
131 // Derivatives of metric coefficients
132 // ----------------------------------
133
134 // ... (d/dX,d/dY,d/dZ)(n_auto) :
135 d_n_auto = n_auto.gradient() ; // (d/dx, d/dy, d/dz)
136 d_n_auto.change_triad(ref_triad) ; // --> (d/dX, d/dY, d/dZ)
137
138 // ... (d/dX,d/dY,d/dZ)(confpsi_auto) :
139 d_confpsi_auto = confpsi_auto.gradient() ; // (d/dx, d/dy, d/dz)
140 d_confpsi_auto.change_triad(ref_triad) ; // --> (d/dX, d/dY, d/dZ)
141
142 // The derived quantities are obsolete
143 // -----------------------------------
144
145 del_deriv() ;
146}
147}
Black hole.
Definition bhole.h:268
const Tenseur & get_n_auto() const
Returns the part of N generated by the hole.
Definition bhole.h:395
const Tenseur & get_psi_auto() const
Returns the part of generated by the hole.
Definition bhole.h:412
const Tenseur & get_shift_auto() const
Returns the part of generated by the hole.
Definition bhole.h:440
Tenseur confpsi_auto
Part of the conformal factor $\Psi$ generated principaly by the star.
Tenseur confpsi_comp
Part of the conformal factor $\Psi$ generated principaly by the companion star.
void update_metric(const Bhole &comp)
Computes metric coefficients from known potentials, when the companion is a black hole.
Tenseur d_n_auto
Gradient of {\tt n_auto} (Cartesian components with respect to {\tt ref_triad}).
Definition et_bin_nsbh.h:93
Tenseur n_auto
Part of the lapse {\it N} generated principaly by the star.
Definition et_bin_nsbh.h:85
Tenseur n_comp
Part of the lapse {\it N} generated principaly by the companion star.
Definition et_bin_nsbh.h:88
Tenseur confpsi
Total conformal factor $\Psi$.
Tenseur d_confpsi_auto
Gradient of {\tt confpsi_auto} (Cartesian components with respect to {\tt ref_triad}...
Tenseur shift_comp
Part of the shift vector generated principaly by the companion star.
Definition etoile.h:898
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
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 nnn
Total lapse function.
Definition etoile.h:512
Tenseur shift
Total shift vector.
Definition etoile.h:515
Tenseur a_car
Total conformal factor .
Definition etoile.h:518
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
Lorene prototypes.
Definition app_hor.h:67