LORENE
bin_bhns_extr.h
1/*
2 * Definition of Lorene class Bin_bhns_extr
3 *
4 */
5
6/*
7 * Copyright (c) 2004-2005 Keisuke Taniguchi
8 *
9 * This file is part of LORENE.
10 *
11 * LORENE is free software; you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License version 2
13 * as published by the Free Software Foundation.
14 *
15 * LORENE is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with LORENE; if not, write to the Free Software
22 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23 *
24 */
25
26#ifndef __BIN_BHNS_EXTR_H_
27#define __BIN_BHNS_EXTR_H_
28
29/*
30 * $Id: bin_bhns_extr.h,v 1.3 2014/10/13 08:52:32 j_novak Exp $
31 * $Log: bin_bhns_extr.h,v $
32 * Revision 1.3 2014/10/13 08:52:32 j_novak
33 * Lorene classes and functions now belong to the namespace Lorene.
34 *
35 * Revision 1.2 2005/02/28 23:02:51 k_taniguchi
36 * Addition of some codes for the conformally flat case
37 *
38 * Revision 1.1 2004/11/30 20:36:00 k_taniguchi
39 * *** empty log message ***
40 *
41 *
42 * $Header: /cvsroot/Lorene/C++/Include/bin_bhns_extr.h,v 1.3 2014/10/13 08:52:32 j_novak Exp $
43 *
44 */
45
46// Lorene header
47#include "et_bin_bhns_extr.h"
48#include "etoile.h"
49
50namespace Lorene {
58
59 // Data :
60 // -----
61 private:
64
67
71 double omega ;
72
74 double separ ;
75
77 double mass_bh ;
78
79 // Derived data :
80 // ------------
81 private :
86 mutable double* p_xa_barycenter_extr ;
87
91 mutable double* p_ya_barycenter_extr ;
92
96 mutable double* p_mass_b_extr ;
97
98 // Constructors - Destructor
99 // -------------------------
100 public:
116 Bin_bhns_extr(Map& mp, int nzet, const Eos& eos,
117 bool irrot, bool relat, bool kerrs, bool multi) ;
118
119 Bin_bhns_extr(const Bin_bhns_extr& ) ;
120
127 Bin_bhns_extr(Map& mp, const Eos& eos, FILE* fich) ;
128
129 ~Bin_bhns_extr() ;
130
131
132 // Memory management
133 // -----------------
134 private:
136 void del_deriv() const ;
137
139 void set_der_0x0() const ;
140
141
142 // Mutators / assignment
143 // ---------------------
144 public:
146 void operator=(const Bin_bhns_extr&) ;
147
150 { del_deriv() ;
151 return star ; } ;
152
154 double& set_omega() { return omega ; } ;
155
157 double& set_separ() { return separ ; } ;
158
160 double& set_mass_bh() { return mass_bh ; } ;
161
162 // Accessors
163 // ---------
164 public:
167 { return star ; } ;
168
170 double get_omega() const { return omega ; } ;
171
175 double get_separ() const { return separ ; } ;
176
178 double get_mass_bh() const { return mass_bh ; } ;
179
180 // Outputs
181 // -------
182 public:
183 void sauve(FILE* ) const ;
184
186 friend ostream& operator<<(ostream& , const Bin_bhns_extr& ) ;
187
189 void display_poly(ostream& ) const ;
190
191 private:
193 ostream& operator>>(ostream& ) const ;
194
195 // Computational routines
196 // ----------------------
197 public:
202 double xa_barycenter_extr() const ;
203
204 // Absolute coordinate Y of the barycenter of the baryon density
206 double ya_barycenter_extr() const ;
207
211 double mass_b_extr() const ;
212
233 void orbit_omega_ks(double fact_omeg_min, double fact_omeg_max) ;
234
255 void orbit_omega_cf(double fact_omeg_min, double fact_omeg_max) ;
256
260 void analytical_omega() ;
261
266 void analytical_shift() ;
267
268};
269ostream& operator<<(ostream& , const Bin_bhns_extr& ) ;
270
271}
272#endif
Cartesian vectorial bases (triads).
Definition base_vect.h:201
Class for computing a Black hole - Neutron star binary system with an extreme mass ratio.
void orbit_omega_ks(double fact_omeg_min, double fact_omeg_max)
Computes the orbital angular velocity {\tt omega} in the Kerr-Schild background metric.
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
double * p_ya_barycenter_extr
Absolute coordinate Y of the barycenter of the baryon density in the Kerr-Schild background metric.
const Et_bin_bhns_extr & get_ns() const
Returns a reference to the neutron star.
double get_omega() const
Returns the orbital angular velocity [{\tt f_unit}].
double get_separ() const
Returns the coordinate separation of the binary system [{\tt r_unit}].
Et_bin_bhns_extr & set_ns()
Read/write of the neutron star.
double * p_mass_b_extr
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat one.
double ya_barycenter_extr() const
in the Kerr-Schild background metric
friend ostream & operator<<(ostream &, const Bin_bhns_extr &)
Display.
double mass_b_extr() const
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat.
Et_bin_bhns_extr star
Neutron star.
void operator=(const Bin_bhns_extr &)
Assignment to another Bin_bhns_extr.
double & set_separ()
Sets the orbital separation [{\tt r_unit}].
double & set_omega()
Sets the orbital angular velocity [{\tt f_unit}].
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
double separ
Absolute orbital separation between two centers of BH and NS.
void analytical_omega()
Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian cas...
void orbit_omega_cf(double fact_omeg_min, double fact_omeg_max)
Computes the orbital angular velocity {\tt omega} in the conformally flat background metric.
double mass_bh
Gravitational mass of BH.
double xa_barycenter_extr() const
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or...
void display_poly(ostream &) const
Display in polytropic units.
void analytical_shift()
Sets some analytical template for the shift vector (via the members {\tt w_shift} and {\tt khi_shift}...
void sauve(FILE *) const
Save in a file.
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<)
double & set_mass_bh()
Sets the gravitational mass of BH [{\tt m_unit}].
~Bin_bhns_extr()
Destructor.
double omega
Angular velocity with respect to an asymptotically inertial observer.
double get_mass_bh() const
Returns the gravitational mass of BH [{\tt m_unit}].
double * p_xa_barycenter_extr
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or...
void del_deriv() const
Deletes all the derived quantities.
Equation of state base class.
Definition eos.h:190
Class for a neutron star in black hole - neutron star binary systems.
Base class for coordinate mappings.
Definition map.h:670
Lorene prototypes.
Definition app_hor.h:64