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 
62 namespace 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 ) {
78  }
79  else {
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 
93  }
94  else {
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 
109  }
110 
111  // Conformal factor
112  // ----------------
113 
114  if ( (hole.get_confo_auto()).get_etat() == ETATZERO ) {
116  }
117  else {
119  confo_comp.import( hole.get_confo_auto() ) ;
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 
144 
147 
150 
151  psi4 = pow(confo_tot, 4.) ;
153 
156 
159 
160  //---------------------------------
161  // Derivative of metric quantities
162  //---------------------------------
163 
165  for (int i=1; i<=3; i++)
167 
169 
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 
178 
180  for (int i=1; i<=3; i++)
181  d_confo_auto.set(i) = confo_auto.deriv(i) ;
182 
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 }
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 shift_tot
Total shift vector.
Definition: star_bhns.h:144
Scalar lapconf_tot
Total lapconf function.
Definition: star_bhns.h:119
Class for stars in black hole-neutron star binaries.
Definition: star_bhns.h:59
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition: tensor.C:490
Map & mp
Mapping associated with the star.
Definition: star.h:180
const Scalar & get_confo_auto() const
Returns the part of the conformal factor generated by the black hole.
Definition: hole_bhns.h:439
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star_bhns.C:344
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
Lorene prototypes.
Definition: app_hor.h:67
Base class for coordinate mappings.
Definition: map.h:688
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Vector shift_auto
Shift vector generated by the star.
Definition: star_bhns.h:138
Tensor d_shift_auto
Derivative of the shift vector generated by the star .
Definition: star_bhns.h:149
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
Tensor field of valence 1.
Definition: vector.h:188
Scalar lapse_tot
Total lapse function.
Definition: star_bhns.h:125
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: vector.C:322
Scalar confo_tot
Total conformal factor.
Definition: star_bhns.h:163
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tensor.h:879
const Scalar & get_lapconf_auto() const
Returns the part of the lapconf function generated by the black hole.
Definition: hole_bhns.h:370
const Map & get_mp() const
Returns the mapping.
Definition: blackhole.h:213
Scalar lapconf_comp
Lapconf function generated by the companion black hole.
Definition: star_bhns.h:116
const Scalar & deriv(int i) const
Returns of *this , where .
Definition: scalar_deriv.C:359
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Vector shift_comp
Shift vector generated by the companion black hole.
Definition: star_bhns.h:141
const Vector & get_shift_auto() const
Returns the part of the shift vector generated by the black hole.
Definition: hole_bhns.h:408
Class for black holes in black hole-neutron star binaries.
Definition: hole_bhns.h:65
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...
void import(const Scalar &ci)
Assignment to another Scalar defined on a different mapping.
Definition: scalar_import.C:71
Scalar lapconf_auto
Lapconf function generated by the star.
Definition: star_bhns.h:113
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:809
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:663
Vector d_confo_auto
Derivative of the conformal factor generated by the star .
Definition: star_bhns.h:168
Scalar confo_auto
Conformal factor generated by the star.
Definition: star_bhns.h:157
Scalar confo_comp
Conformal factor generated by the companion black hole.
Definition: star_bhns.h:160
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tensor.C:528
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:506
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
Scalar lapse_auto
Lapse function generated by the "star".
Definition: star_bhns.h:122
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: tensor.C:935