LORENE
bin_bhns_rotaxis.C
1 /*
2  * Methods of class Bin_bhns to compute the location of the rotation axis
3  *
4  * (see file bin_bhns.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2005-2006 Keisuke Taniguchi
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License version 2
15  * as published by the Free Software Foundation.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 
30 /*
31  * $Id: bin_bhns_rotaxis.C,v 1.4 2016/12/05 16:17:45 j_novak Exp $
32  * $Log: bin_bhns_rotaxis.C,v $
33  * Revision 1.4 2016/12/05 16:17:45 j_novak
34  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35  *
36  * Revision 1.3 2014/10/13 08:52:41 j_novak
37  * Lorene classes and functions now belong to the namespace Lorene.
38  *
39  * Revision 1.2 2014/10/06 15:13:00 j_novak
40  * Modified #include directives to use c++ syntax.
41  *
42  * Revision 1.1 2007/06/22 01:10:47 k_taniguchi
43  * *** empty log message ***
44  *
45  *
46  * $Header: /cvsroot/Lorene/C++/Source/Bin_bhns/bin_bhns_rotaxis.C,v 1.4 2016/12/05 16:17:45 j_novak Exp $
47  *
48  */
49 
50 // C++ headers
51 //#include <>
52 
53 // C headers
54 #include <cmath>
55 
56 // Lorene headers
57 #include "bin_bhns.h"
58 #include "unites.h"
59 
60 namespace Lorene {
61 void Bin_bhns::rotation_axis_x(double rot_exp_x) {
62 
63  using namespace Unites ;
64 
65  double momunit = (hole.get_mass_bh()+star.mass_g_bhns())*omega*separ ;
66 
67  double error_y = line_mom_bhns()(1) / momunit ;
68 
69  if (error_y >= 1.) {
70  cout << "Bin_bhns::rotation_axis:" << endl ;
71  cout << " !!! WARNING : error_y is larger than +1 !!!" << endl ;
72  error_y *= 0.1 ;
73  }
74 
75  // Sets X_BH and X_NS
76  // ------------------
77 
78  double gg = pow( (2.-error_y)/(2.-2.*error_y), rot_exp_x) ;
79 
80  double xbh_old = hole.get_mp().get_ori_x() ;
81 
82  cout << "Bin_bhns::rotation_axis:" << endl ;
83  cout << " error_y : " << error_y << " gg : " << gg << endl ;
84  cout << " old X_BH : " << hole.get_mp().get_ori_x() / km << " [km]"
85  << " old X_NS : " << star.get_mp().get_ori_x() / km << " [km]"
86  << endl ;
87 
88  double xbh_new = xbh_old * gg ;
89  double xns_new = xbh_new + separ ;
90 
91  cout << " new X_BH : " << xbh_new / km << " [km]"
92  << " new X_NS : " << xns_new / km << " [km]"
93  << endl ;
94 
95  double yns_old = star.get_mp().get_ori_y() ;
96 
97  (hole.set_mp()).set_ori(xbh_new, 0., 0.) ;
98  (star.set_mp()).set_ori(xns_new, yns_old, 0.) ;
99 
100  set_x_rot() = 0. ;
101 
102 }
103 
104 void Bin_bhns::rotation_axis_y(double thres_rot, double rot_exp_y,
105  double fact) {
106 
107  using namespace Unites ;
108 
109  double momunit = (hole.get_mass_bh()+star.mass_g_bhns())*omega*separ ;
110 
111  double error_x = line_mom_bhns()(0) / momunit ;
112 
113  if (error_x <= -1.) {
114  cout << "Bin_bhns::rotation_axis:" << endl ;
115  cout << " !!! WARNING : error_x is smaller than -1 !!!" << endl ;
116  error_x *= 0.1 ;
117  }
118 
119  // Sets Y_NS
120  // ---------
121 
122  double ff = pow( (2.+error_x)/(2.+2.*error_x), rot_exp_y) ;
123 
124  double yns_old = star.get_mp().get_ori_y() ;
125 
126  if ( fabs(error_x) < thres_rot ) {
127  cout << "Bin_bhns::rotation_axis:" << endl ;
128  cout << " ff is set to 1 because error_x is smaller than" << endl ;
129  cout << " the threshold value (" << thres_rot << ")" << endl ;
130 
131  ff = 1. ;
132  }
133 
134  cout << "Local center of mass of NS:" << endl ;
135  cout << " X_CM : " << xa_barycenter() / km << " [km]"
136  << " Y_CM : " << ya_barycenter() / km << " [km]" << endl ;
137 
138  cout << "Bin_bhns::rotation_axis:" << endl ;
139  cout << " error_x : " << error_x << " ff : " << ff << endl ;
140  cout << " old Y_BH : " << hole.get_mp().get_ori_y() / km << " [km]"
141  << " old Y_NS : " << star.get_mp().get_ori_y() / km << " [km]"
142  << endl ;
143 
144  double aa = fact * separ ;
145  double yns_new = yns_old + aa * (1. - ff) ;
146 
147  cout << " new Y_BH : " << 0. / km << " [km]"
148  << " new Y_NS : " << yns_new / km << " [km]"
149  << endl ;
150 
151  double xbh_old = hole.get_mp().get_ori_x() ;
152  double xns_old = star.get_mp().get_ori_x() ;
153 
154  (hole.set_mp()).set_ori(xbh_old, 0., 0.) ;
155  (star.set_mp()).set_ori(xns_old, yns_new, 0.) ;
156 
157  set_y_rot() = 0. ;
158 
159 }
160 }
double get_mass_bh() const
Returns the gravitational mass of BH [{ m_unit}].
Definition: blackhole.h:221
Hole_bhns hole
Black hole.
Definition: bin_bhns.h:72
const Tbl & line_mom_bhns() const
Total linear momentum.
double ya_barycenter() const
Absolute coordinate Y of the barycenter of the baryon density.
double & set_x_rot()
Sets the absolute coordinate X of the rotation axis [{ r_unit}].
Definition: bin_bhns.h:215
double get_ori_y() const
Returns the y coordinate of the origin.
Definition: map.h:782
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
double xa_barycenter() const
Absolute coordinate X of the barycenter of the baryon density.
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:780
const Map & get_mp() const
Returns the mapping.
Definition: star.h:355
double separ
Absolute orbital separation between two centers of BH and NS.
Definition: bin_bhns.h:83
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition: bin_bhns.h:80
Map & set_mp()
Read/write of the mapping.
Definition: blackhole.h:204
Star_bhns star
Neutron star.
Definition: bin_bhns.h:75
const Map & get_mp() const
Returns the mapping.
Definition: blackhole.h:213
void rotation_axis_x(double rot_exp_x)
Computes the position of the rotation axis X.
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
double & set_y_rot()
Sets the absolute coordinate Y of the rotation axis [{ r_unit}].
Definition: bin_bhns.h:218
void rotation_axis_y(double thres_rot, double rot_exp_y, double fact)
Computes the position of the rotation axis Y.
Map & set_mp()
Read/write of the mapping.
Definition: star.h:322