LORENE
hole_bhns_upmetr_der.C
1
/*
2
* Method of class Hole_bhns to compute the derivative of metric quantities
3
* from the companion neutron-star
4
*
5
* (see file hole_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
char
hole_bhns_upmetr_der_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_upmetr_der.C,v 1.3 2014/10/13 08:53:00 j_novak Exp $"
;
30
31
/*
32
* $Id: hole_bhns_upmetr_der.C,v 1.3 2014/10/13 08:53:00 j_novak Exp $
33
* $Log: hole_bhns_upmetr_der.C,v $
34
* Revision 1.3 2014/10/13 08:53:00 j_novak
35
* Lorene classes and functions now belong to the namespace Lorene.
36
*
37
* Revision 1.2 2008/05/15 19:09:02 k_taniguchi
38
* Change of some parameters.
39
*
40
* Revision 1.1 2007/06/22 01:25:50 k_taniguchi
41
* *** empty log message ***
42
*
43
*
44
* $Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_upmetr_der.C,v 1.3 2014/10/13 08:53:00 j_novak Exp $
45
*
46
*/
47
48
// C++ headers
49
//#include <>
50
51
// C headers
52
//#include <>
53
54
// Lorene headers
55
#include "hole_bhns.h"
56
#include "star_bhns.h"
57
#include "utilitaires.h"
58
59
namespace
Lorene
{
60
void
Hole_bhns::update_met_der_comp_bhns
(
const
Star_bhns
& star) {
61
62
// Computation of d_lapconf_comp
63
// -----------------------------
64
65
if
( (star.
get_d_lapconf_auto
())(1).get_etat() == ETATZERO ) {
66
assert( (star.
get_d_lapconf_auto
())(2).get_etat() == ETATZERO ) ;
67
assert( (star.
get_d_lapconf_auto
())(3).get_etat() == ETATZERO ) ;
68
69
d_lapconf_comp
.set_etat_zero() ;
70
}
71
else
{
72
d_lapconf_comp
.set_etat_qcq() ;
73
Vector
comp_dlapconf( star.
get_d_lapconf_auto
() ) ;
74
comp_dlapconf.
dec_dzpuis
(2) ;
// dzpuis : 2 -> 0 for import
75
76
(
d_lapconf_comp
.set(1)).import( comp_dlapconf(1) ) ;
77
(
d_lapconf_comp
.set(2)).import( comp_dlapconf(2) ) ;
78
(
d_lapconf_comp
.set(3)).import( comp_dlapconf(3) ) ;
79
80
d_lapconf_comp
.std_spectral_base() ;
81
d_lapconf_comp
.inc_dzpuis(2) ;
// dzpuis : 0 -> 2
82
}
83
84
85
// Computation of d_shift_comp
86
// ---------------------------
87
88
if
( (star.
get_d_shift_auto
())(1,2).get_etat() == ETATZERO ) {
89
assert( (star.
get_d_shift_auto
())(1,1).get_etat() == ETATZERO ) ;
90
assert( (star.
get_d_shift_auto
())(1,3).get_etat() == ETATZERO ) ;
91
92
d_shift_comp
.set_etat_zero() ;
93
}
94
else
{
95
96
d_shift_comp
.set_etat_qcq() ;
97
Tensor
comp_dshift( star.
get_d_shift_auto
() ) ;
98
comp_dshift.
dec_dzpuis
(2) ;
// dzpuis : 2 -> 0 for import
99
100
(
d_shift_comp
.set(1,1)).import( comp_dshift(1,1) ) ;
101
(
d_shift_comp
.set(1,2)).import( comp_dshift(1,2) ) ;
102
(
d_shift_comp
.set(1,3)).import( comp_dshift(1,3) ) ;
103
(
d_shift_comp
.set(2,1)).import( comp_dshift(2,1) ) ;
104
(
d_shift_comp
.set(2,2)).import( comp_dshift(2,2) ) ;
105
(
d_shift_comp
.set(2,3)).import( comp_dshift(2,3) ) ;
106
(
d_shift_comp
.set(3,1)).import( comp_dshift(3,1) ) ;
107
(
d_shift_comp
.set(3,2)).import( comp_dshift(3,2) ) ;
108
(
d_shift_comp
.set(3,3)).import( comp_dshift(3,3) ) ;
109
110
d_shift_comp
.std_spectral_base() ;
111
d_shift_comp
.inc_dzpuis(2) ;
// dzpuis : 0 -> 2
112
}
113
114
115
// Computation of d_confo_comp
116
// ---------------------------
117
118
if
( (star.
get_d_confo_auto
())(1).get_etat() == ETATZERO ) {
119
assert( (star.
get_d_confo_auto
())(2).get_etat() == ETATZERO ) ;
120
assert( (star.
get_d_confo_auto
())(3).get_etat() == ETATZERO ) ;
121
122
d_confo_comp
.set_etat_zero() ;
123
}
124
else
{
125
d_confo_comp
.set_etat_qcq() ;
126
Vector
comp_dconfo( star.
get_d_confo_auto
() ) ;
127
comp_dconfo.
dec_dzpuis
(2) ;
// dzpuis : 2 -> 0 for import
128
129
(
d_confo_comp
.set(1)).import( comp_dconfo(1) ) ;
130
(
d_confo_comp
.set(2)).import( comp_dconfo(2) ) ;
131
(
d_confo_comp
.set(3)).import( comp_dconfo(3) ) ;
132
133
d_confo_comp
.std_spectral_base() ;
134
d_confo_comp
.inc_dzpuis(2) ;
// dzpuis : 0 -> 2
135
}
136
137
138
// The derived quantities are obsolete
139
// -----------------------------------
140
141
del_deriv
() ;
142
143
}
144
}
Lorene::Hole_bhns::del_deriv
virtual void del_deriv() const
Deletes all the derived quantities.
Definition
hole_bhns.C:392
Lorene::Hole_bhns::d_lapconf_comp
Vector d_lapconf_comp
Derivative of the lapconf function generated by the companion star.
Definition
hole_bhns.h:123
Lorene::Hole_bhns::update_met_der_comp_bhns
void update_met_der_comp_bhns(const Star_bhns &star)
Computes derivative of metric quantities from the companion neutron star.
Definition
hole_bhns_upmetr_der.C:60
Lorene::Hole_bhns::d_confo_comp
Vector d_confo_comp
Derivative of the conformal factor generated by the companion star.
Definition
hole_bhns.h:185
Lorene::Hole_bhns::d_shift_comp
Tensor d_shift_comp
Derivative of the shift vector generated by the companion star.
Definition
hole_bhns.h:154
Lorene::Star_bhns
Class for stars in black hole-neutron star binaries.
Definition
star_bhns.h:59
Lorene::Star_bhns::get_d_lapconf_auto
const Vector & get_d_lapconf_auto() const
Returns the derivative of the lapse function generated by the star.
Definition
star_bhns.h:356
Lorene::Star_bhns::get_d_shift_auto
const Tensor & get_d_shift_auto() const
Returns the derivative of the shift vector generated by the star.
Definition
star_bhns.h:375
Lorene::Star_bhns::get_d_confo_auto
const Vector & get_d_confo_auto() const
Returns the derivative of the conformal factor generated by the star.
Definition
star_bhns.h:396
Lorene::Tensor
Tensor handling.
Definition
tensor.h:288
Lorene::Vector
Tensor field of valence 1.
Definition
vector.h:188
Lorene::Tensor::dec_dzpuis
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition
tensor.C:808
Lorene
Lorene prototypes.
Definition
app_hor.h:64
C++
Source
Hole_bhns
hole_bhns_upmetr_der.C
Generated by
1.15.0