LORENE
star_bhns_upmetr_der.C
1/*
2 * Method of class Star_bhns compute the derivative of metric quantities
3 * from the companion black-hole
4 *
5 * (see file star_bhns.h for documentation).
6 *
7 */
8
9/*
10 * Copyright (c) 2006-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
29char star_bhns_upmetr_der_C[] = "$Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_upmetr_der.C,v 1.3 2014/10/13 08:53:41 j_novak Exp $" ;
30
31/*
32 * $Id: star_bhns_upmetr_der.C,v 1.3 2014/10/13 08:53:41 j_novak Exp $
33 * $Log: star_bhns_upmetr_der.C,v $
34 * Revision 1.3 2014/10/13 08:53:41 j_novak
35 * Lorene classes and functions now belong to the namespace Lorene.
36 *
37 * Revision 1.2 2008/05/15 19:18:17 k_taniguchi
38 * Change of some parameters.
39 *
40 * Revision 1.1 2007/06/22 01:32:55 k_taniguchi
41 * *** empty log message ***
42 *
43 *
44 * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_upmetr_der.C,v 1.3 2014/10/13 08:53:41 j_novak Exp $
45 *
46 */
47
48// C++ headers
49//#include <>
50
51// C headers
52//#include <>
53
54// Lorene headers
55#include "star_bhns.h"
56#include "hole_bhns.h"
57#include "utilitaires.h"
58
59namespace Lorene {
61
62 // Computation of d_lapconf_comp
63 // -----------------------------
64
65 if ( (hole.get_d_lapconf_auto())(1).get_etat() == ETATZERO ) {
66 assert( (hole.get_d_lapconf_auto())(2).get_etat() == ETATZERO ) ;
67 assert( (hole.get_d_lapconf_auto())(3).get_etat() == ETATZERO ) ;
68
70 }
71 else {
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
81 d_lapconf_comp.inc_dzpuis(2) ; // dzpuis : 0 -> 2
82 }
83
84
85 // Computation of d_shift_comp
86 // ---------------------------
87
88 if ( (hole.get_d_shift_auto())(1,2).get_etat() == ETATZERO ) {
89 assert( (hole.get_d_shift_auto())(1,1).get_etat() == ETATZERO ) ;
90 assert( (hole.get_d_shift_auto())(1,3).get_etat() == ETATZERO ) ;
91
93 }
94 else {
95
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
111 d_shift_comp.inc_dzpuis(2) ; // dzpuis : 0 -> 2
112 }
113
114
115 // Computation of d_confo_comp
116 // ---------------------------
117
118 if ( (hole.get_d_confo_auto())(1).get_etat() == ETATZERO ) {
119 assert( (hole.get_d_confo_auto())(2).get_etat() == ETATZERO ) ;
120 assert( (hole.get_d_confo_auto())(3).get_etat() == ETATZERO ) ;
121
123 }
124 else {
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
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}
Time evolution with partial storage (*** under development ***).
Definition evolution.h:371
Class for black holes in black hole-neutron star binaries.
Definition hole_bhns.h:65
const Vector & get_d_confo_auto() const
Returns the derivative of the conformal factor generated by the black hole.
Definition hole_bhns.h:457
const Tensor & get_d_shift_auto() const
Returns the derivative of the shift vector generated by the black hole.
Definition hole_bhns.h:426
const Vector & get_d_lapconf_auto() const
Returns the derivative of the lapconf function generated by the black hole.
Definition hole_bhns.h:397
virtual void del_deriv() const
Deletes all the derived quantities.
Definition star_bhns.C:341
Vector d_confo_comp
Derivative of the conformal factor generated by the companion black hole.
Definition star_bhns.h:173
Vector d_lapconf_comp
Derivative of the lapconf function generated by the companion black hole.
Definition star_bhns.h:135
void update_met_der_comp_bhns(const Hole_bhns &hole)
Computes derivative of metric quantities from the companion black hole.
Tensor d_shift_comp
Derivative of the shift vector generated by the companion black hole.
Definition star_bhns.h:154
Tensor handling.
Definition tensor.h:288
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:316
Scalar & set(int)
Read/write access to a component.
Definition vector.C:296
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:816
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition tensor.C:481
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition tensor.C:497
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition tensor.C:654
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition tensor.C:926
Lorene prototypes.
Definition app_hor.h:64