LORENE
star_bhns_upmetr.C
1/*
2 * Method of class Star_bhns to compute metric quantities from
3 * the companion black-hole and total metric quantities
4 *
5 * (see file star_bhns.h for documentation).
6 *
7 */
8
9/*
10 * Copyright (c) 2005,2007 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: star_bhns_upmetr.C,v 1.4 2016/12/05 16:18:16 j_novak Exp $
33 * $Log: star_bhns_upmetr.C,v $
34 * Revision 1.4 2016/12/05 16:18:16 j_novak
35 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36 *
37 * Revision 1.3 2014/10/13 08:53:41 j_novak
38 * Lorene classes and functions now belong to the namespace Lorene.
39 *
40 * Revision 1.2 2008/05/15 19:17:39 k_taniguchi
41 * Change of some parameters.
42 *
43 * Revision 1.1 2007/06/22 01:32:37 k_taniguchi
44 * *** empty log message ***
45 *
46 *
47 * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_upmetr.C,v 1.4 2016/12/05 16:18:16 j_novak Exp $
48 *
49 */
50
51// C++ headers
52//#include <>
53
54// C headers
55//#include <>
56
57// Lorene headers
58#include "star_bhns.h"
59#include "hole_bhns.h"
60#include "utilitaires.h"
61
62namespace Lorene {
64 const Star_bhns& star_prev,
65 double relax_met_comp) {
66
67 //-----------------------------------------------------
68 // Computation of quantities coming from the companion
69 //-----------------------------------------------------
70
71 const Map& mp_bh (hole.get_mp()) ;
72
73 // Lapconf function
74 // ----------------
75
76 if ( (hole.get_lapconf_auto()).get_etat() == ETATZERO ) {
77 lapconf_comp.set_etat_zero() ;
78 }
79 else {
80 lapconf_comp.set_etat_qcq() ;
81 lapconf_comp.import( hole.get_lapconf_auto() ) ;
82 lapconf_comp.std_spectral_base() ;
83 }
84
85 // Shift vector
86 // ------------
87
88 if ( (hole.get_shift_auto())(2).get_etat() == ETATZERO ) {
89 assert( (hole.get_shift_auto())(1).get_etat() == ETATZERO ) ;
90 assert( (hole.get_shift_auto())(3).get_etat() == ETATZERO ) ;
91
92 shift_comp.set_etat_zero() ;
93 }
94 else {
95 shift_comp.set_etat_qcq() ;
96 shift_comp.set_triad(mp.get_bvect_cart()) ;
97
98 Vector comp_shift(hole.get_shift_auto()) ;
99 comp_shift.change_triad(mp_bh.get_bvect_cart()) ;
100 comp_shift.change_triad(mp.get_bvect_cart()) ;
101
102 assert( *(shift_comp.get_triad()) == *(comp_shift.get_triad()) ) ;
103
104 (shift_comp.set(1)).import( comp_shift(1) ) ;
105 (shift_comp.set(2)).import( comp_shift(2) ) ;
106 (shift_comp.set(3)).import( comp_shift(3) ) ;
107
108 shift_comp.std_spectral_base() ;
109 }
110
111 // Conformal factor
112 // ----------------
113
114 if ( (hole.get_confo_auto()).get_etat() == ETATZERO ) {
115 confo_comp.set_etat_zero() ;
116 }
117 else {
118 confo_comp.set_etat_qcq() ;
119 confo_comp.import( hole.get_confo_auto() ) ;
120 confo_comp.std_spectral_base() ;
121 }
122
123 //----------------------------------------------------
124 // Relaxation on lapconf_comp, shift_comp, confo_comp
125 //----------------------------------------------------
126
127 double relax_met_comp_jm1 = 1. - relax_met_comp ;
128
129 lapconf_comp = relax_met_comp * lapconf_comp
130 + relax_met_comp_jm1 * (star_prev.lapconf_comp) ;
131
132 shift_comp = relax_met_comp * shift_comp
133 + relax_met_comp_jm1 * (star_prev.shift_comp) ;
134
135 confo_comp = relax_met_comp * confo_comp
136 + relax_met_comp_jm1 * (star_prev.confo_comp) ;
137
138 //-------------------------
139 // Total metric quantities
140 //-------------------------
141
143 lapconf_tot.std_spectral_base() ;
144
146 shift_tot.std_spectral_base() ;
147
149 confo_tot.std_spectral_base() ;
150
151 psi4 = pow(confo_tot, 4.) ;
152 psi4.std_spectral_base() ;
153
155 lapse_auto.std_spectral_base() ;
156
158 lapse_tot.std_spectral_base() ;
159
160 //---------------------------------
161 // Derivative of metric quantities
162 //---------------------------------
163
164 d_lapconf_auto.set_etat_qcq() ;
165 for (int i=1; i<=3; i++)
166 d_lapconf_auto.set(i) = lapconf_auto.deriv(i) ;
167
168 d_lapconf_auto.std_spectral_base() ;
169
170 d_shift_auto.set_etat_qcq() ;
171 for (int i=1; i<=3; i++) {
172 for (int j=1; j<=3; j++) {
173 d_shift_auto.set(i,j) = shift_auto(j).deriv(i) ;
174 }
175 }
176
177 d_shift_auto.std_spectral_base() ;
178
179 d_confo_auto.set_etat_qcq() ;
180 for (int i=1; i<=3; i++)
181 d_confo_auto.set(i) = confo_auto.deriv(i) ;
182
183 d_confo_auto.std_spectral_base() ;
184
185 // ... extrinsic curvature (taij_auto and taij_quad_auto)
186 // extr_curv_bhns() ;
187
188 // The derived quantities are obsolete
189 // -----------------------------------
190
191 del_deriv() ;
192
193}
194}
const Map & get_mp() const
Returns the mapping.
Definition blackhole.h:213
Class for black holes in black hole-neutron star binaries.
Definition hole_bhns.h:65
const Scalar & get_lapconf_auto() const
Returns the part of the lapconf function generated by the black hole.
Definition hole_bhns.h:370
const Scalar & get_confo_auto() const
Returns the part of the conformal factor generated by the black hole.
Definition hole_bhns.h:439
const Vector & get_shift_auto() const
Returns the part of the shift vector generated by the black hole.
Definition hole_bhns.h:408
Vector shift_tot
Total shift vector.
Definition star_bhns.h:144
Scalar lapconf_auto
Lapconf function generated by the star.
Definition star_bhns.h:113
Scalar confo_comp
Conformal factor generated by the companion black hole.
Definition star_bhns.h:160
Tensor d_shift_auto
Derivative of the shift vector generated by the star .
Definition star_bhns.h:149
Scalar lapse_tot
Total lapse function.
Definition star_bhns.h:125
virtual void del_deriv() const
Deletes all the derived quantities.
Definition star_bhns.C:344
Scalar confo_tot
Total conformal factor.
Definition star_bhns.h:163
Scalar psi4
Fourth power of the total conformal factor.
Definition star_bhns.h:176
Vector d_lapconf_auto
Derivative of the lapconf function generated by the star .
Definition star_bhns.h:130
Vector d_confo_auto
Derivative of the conformal factor generated by the star .
Definition star_bhns.h:168
Vector shift_comp
Shift vector generated by the companion black hole.
Definition star_bhns.h:141
Scalar lapconf_tot
Total lapconf function.
Definition star_bhns.h:119
Star_bhns(Map &mp_i, int nzet_i, const Eos &eos_i, bool irrot_i)
Standard constructor.
Definition star_bhns.C:76
Scalar lapse_auto
Lapse function generated by the "star".
Definition star_bhns.h:122
Scalar lapconf_comp
Lapconf function generated by the companion black hole.
Definition star_bhns.h:116
void update_metric_bhns(const Hole_bhns &hole, const Star_bhns &star_prev, double relax)
Computes metric coefficients from known potentials with relaxation when the companion is a black hole...
Vector shift_auto
Shift vector generated by the star.
Definition star_bhns.h:138
Scalar confo_auto
Conformal factor generated by the star.
Definition star_bhns.h:157
Map & mp
Mapping associated with the star.
Definition star.h:180
Tensor field of valence 1.
Definition vector.h:188
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition tensor.h:879
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142