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 
30 
31 /*
32  * $Id: hole_bhns_upmetr_der.C,v 1.4 2016/12/05 16:17:55 j_novak Exp $
33  * $Log: hole_bhns_upmetr_der.C,v $
34  * Revision 1.4 2016/12/05 16:17:55 j_novak
35  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36  *
37  * Revision 1.3 2014/10/13 08:53:00 j_novak
38  * Lorene classes and functions now belong to the namespace Lorene.
39  *
40  * Revision 1.2 2008/05/15 19:09:02 k_taniguchi
41  * Change of some parameters.
42  *
43  * Revision 1.1 2007/06/22 01:25:50 k_taniguchi
44  * *** empty log message ***
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_upmetr_der.C,v 1.4 2016/12/05 16:17:55 j_novak Exp $
48  *
49  */
50 
51 // C++ headers
52 //#include <>
53 
54 // C headers
55 //#include <>
56 
57 // Lorene headers
58 #include "hole_bhns.h"
59 #include "star_bhns.h"
60 #include "utilitaires.h"
61 
62 namespace Lorene {
64 
65  // Computation of d_lapconf_comp
66  // -----------------------------
67 
68  if ( (star.get_d_lapconf_auto())(1).get_etat() == ETATZERO ) {
69  assert( (star.get_d_lapconf_auto())(2).get_etat() == ETATZERO ) ;
70  assert( (star.get_d_lapconf_auto())(3).get_etat() == ETATZERO ) ;
71 
73  }
74  else {
76  Vector comp_dlapconf( star.get_d_lapconf_auto() ) ;
77  comp_dlapconf.dec_dzpuis(2) ; // dzpuis : 2 -> 0 for import
78 
79  (d_lapconf_comp.set(1)).import( comp_dlapconf(1) ) ;
80  (d_lapconf_comp.set(2)).import( comp_dlapconf(2) ) ;
81  (d_lapconf_comp.set(3)).import( comp_dlapconf(3) ) ;
82 
84  d_lapconf_comp.inc_dzpuis(2) ; // dzpuis : 0 -> 2
85  }
86 
87 
88  // Computation of d_shift_comp
89  // ---------------------------
90 
91  if ( (star.get_d_shift_auto())(1,2).get_etat() == ETATZERO ) {
92  assert( (star.get_d_shift_auto())(1,1).get_etat() == ETATZERO ) ;
93  assert( (star.get_d_shift_auto())(1,3).get_etat() == ETATZERO ) ;
94 
96  }
97  else {
98 
100  Tensor comp_dshift( star.get_d_shift_auto() ) ;
101  comp_dshift.dec_dzpuis(2) ; // dzpuis : 2 -> 0 for import
102 
103  (d_shift_comp.set(1,1)).import( comp_dshift(1,1) ) ;
104  (d_shift_comp.set(1,2)).import( comp_dshift(1,2) ) ;
105  (d_shift_comp.set(1,3)).import( comp_dshift(1,3) ) ;
106  (d_shift_comp.set(2,1)).import( comp_dshift(2,1) ) ;
107  (d_shift_comp.set(2,2)).import( comp_dshift(2,2) ) ;
108  (d_shift_comp.set(2,3)).import( comp_dshift(2,3) ) ;
109  (d_shift_comp.set(3,1)).import( comp_dshift(3,1) ) ;
110  (d_shift_comp.set(3,2)).import( comp_dshift(3,2) ) ;
111  (d_shift_comp.set(3,3)).import( comp_dshift(3,3) ) ;
112 
114  d_shift_comp.inc_dzpuis(2) ; // dzpuis : 0 -> 2
115  }
116 
117 
118  // Computation of d_confo_comp
119  // ---------------------------
120 
121  if ( (star.get_d_confo_auto())(1).get_etat() == ETATZERO ) {
122  assert( (star.get_d_confo_auto())(2).get_etat() == ETATZERO ) ;
123  assert( (star.get_d_confo_auto())(3).get_etat() == ETATZERO ) ;
124 
126  }
127  else {
129  Vector comp_dconfo( star.get_d_confo_auto() ) ;
130  comp_dconfo.dec_dzpuis(2) ; // dzpuis : 2 -> 0 for import
131 
132  (d_confo_comp.set(1)).import( comp_dconfo(1) ) ;
133  (d_confo_comp.set(2)).import( comp_dconfo(2) ) ;
134  (d_confo_comp.set(3)).import( comp_dconfo(3) ) ;
135 
137  d_confo_comp.inc_dzpuis(2) ; // dzpuis : 0 -> 2
138  }
139 
140 
141  // The derived quantities are obsolete
142  // -----------------------------------
143 
144  del_deriv() ;
145 
146 }
147 }
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
Lorene prototypes.
Definition: app_hor.h:67
Tensor d_shift_comp
Derivative of the shift vector generated by the companion star.
Definition: hole_bhns.h:154
const Vector & get_d_confo_auto() const
Returns the derivative of the conformal factor generated by the star.
Definition: star_bhns.h:396
Tensor field of valence 1.
Definition: vector.h:188
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: vector.C:322
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:817
Vector d_lapconf_comp
Derivative of the lapconf function generated by the companion star.
Definition: hole_bhns.h:123
Tensor handling.
Definition: tensor.h:294
const Vector & get_d_lapconf_auto() const
Returns the derivative of the lapse function generated by the star.
Definition: star_bhns.h:356
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:825
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:663
Vector d_confo_comp
Derivative of the conformal factor generated by the companion star.
Definition: hole_bhns.h:185
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
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: tensor.C:935
const Tensor & get_d_shift_auto() const
Returns the derivative of the shift vector generated by the star.
Definition: star_bhns.h:375
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: hole_bhns.C:395
void update_met_der_comp_bhns(const Star_bhns &star)
Computes derivative of metric quantities from the companion neutron star.