Binary parameter class.
More...
#include <binary_tree.h>
|
template<class Tpi , class Tpj > |
void | particleToSemiEccPeriod (const Tpi &_p1, const Tpj &_p2, const Float _G) |
| calcualte semi, ecc and period More...
|
|
Float | calcEccAnomaly (const Float _r) |
| calculate eccentric anomaly from separation (0-pi) More...
|
|
template<class Tptcl > |
void | calcOrbit (const Tptcl &_p1, const Tptcl &_p2, const Float _G) |
| calculate kepler Orbit from particles More...
|
|
template<class Tptcl > |
void | calcParticles (Tptcl &_p1, Tptcl &_p2, const Float &_G) |
| calculate two components from kepler Orbit More...
|
|
void | calcSemiFromPeriod (const Float &_G) |
| from period calculate semi-major axis More...
|
|
template<class Tptcl > |
void | calcParticlesEcca (Tptcl &_p1, Tptcl &_p2, const Float _ecca, const Float _G) const |
| calculate two components from kepler Orbit with input eccentricity anomaly More...
|
|
void | rotateToOriginalFrame (Float *_vec) |
| rotate position vector from binary rest-frame to original frame based on three angles More...
|
|
void | evolve (const Float _dt) |
| Solve kepler motion for dt. More...
|
|
void | print (std::ostream &_os) const |
|
void | printColumn (std::ostream &_fout, const int _width=20) |
| print data of class members using column style More...
|
|
void | writeBinary (FILE *_fp) const |
| write class data to file with binary format More...
|
|
void | readBinary (FILE *_fin) |
| read class data to file with binary format More...
|
|
void | writeAscii (std::ostream &_fout) const |
| write class data to file with ASCII format More...
|
|
void | readAscii (std::istream &_fin) |
| read class data to file with ASCII format More...
|
|
|
template<class Tptcl > |
static void | orbitToParticle (Tptcl &_p1, Tptcl &_p2, const Binary &_bin, const Float &_ecca, const Float _G) |
| Orbit to position and velocity. More...
|
|
template<class Tptcl > |
static void | particleToOrbit (Binary &_bin, const Tptcl &_p1, const Tptcl &_p2, const Float _G) |
| position velocity to orbit More...
|
|
static Float | periodToSemi (const Float &_period, const Float &_mtot, const Float &_G) |
| from period calculate semi-major axis More...
|
|
static Float | semiToPeriod (const Float &_semi, const Float &_mtot, const Float &_G) |
| from semi-major axis calculate period; More...
|
|
template<class Tpi , class Tpj > |
static void | particleToSemiEcc (Float &_semi, Float &_ecc, Float &_r, Float &_rv, const Tpi &_p1, const Tpj &_p2, const Float _G) |
| position velocity to orbit semi-major axis and eccentricity More...
|
|
static Float | calcMeanAnomaly (const Float _ecca, const Float _ecc) |
| calculate mean anomaly from eccentric anomaly (0-pi) More...
|
|
static Float | calcEccAnomaly (const Float _mean_anomaly, const Float _ecc) |
| calculate eccentric anomaly from mean anomaly More...
|
|
static void | solveKepler (Binary &_bin, const Float _dt) |
| solve kepler orbit after dt More...
|
|
static void | printColumnTitle (std::ostream &_fout, const int _width=20) |
| print titles of class members using column style More...
|
|
◆ calcEccAnomaly() [1/2]
static Float COMM::Binary::calcEccAnomaly |
( |
const Float |
_mean_anomaly, |
|
|
const Float |
_ecc |
|
) |
| |
|
inlinestatic |
calculate eccentric anomaly from mean anomaly
◆ calcEccAnomaly() [2/2]
Float COMM::Binary::calcEccAnomaly |
( |
const Float |
_r | ) |
|
|
inline |
calculate eccentric anomaly from separation (0-pi)
- Parameters
-
- Returns
- eccentric anomaly
◆ calcMeanAnomaly()
static Float COMM::Binary::calcMeanAnomaly |
( |
const Float |
_ecca, |
|
|
const Float |
_ecc |
|
) |
| |
|
inlinestatic |
calculate mean anomaly from eccentric anomaly (0-pi)
- Parameters
-
[in] | _ecca | eccentric anomaly |
[in] | _ecc | eccentricity |
- Returns
- mean anomaly
◆ calcOrbit()
template<class Tptcl >
void COMM::Binary::calcOrbit |
( |
const Tptcl & |
_p1, |
|
|
const Tptcl & |
_p2, |
|
|
const Float |
_G |
|
) |
| |
|
inline |
calculate kepler Orbit from particles
- Parameters
-
[out] | _p1 | particle 1 |
[out] | _p2 | particle 2 |
[in] | _G | gravitational constant |
◆ calcParticles()
template<class Tptcl >
void COMM::Binary::calcParticles |
( |
Tptcl & |
_p1, |
|
|
Tptcl & |
_p2, |
|
|
const Float & |
_G |
|
) |
| |
|
inline |
calculate two components from kepler Orbit
- Parameters
-
[out] | _p1 | particle 1 |
[out] | _p2 | particle 2 |
[in] | _G | gravitational constant |
◆ calcParticlesEcca()
template<class Tptcl >
void COMM::Binary::calcParticlesEcca |
( |
Tptcl & |
_p1, |
|
|
Tptcl & |
_p2, |
|
|
const Float |
_ecca, |
|
|
const Float |
_G |
|
) |
| const |
|
inline |
calculate two components from kepler Orbit with input eccentricity anomaly
- Parameters
-
[out] | _p1 | particle 1 |
[out] | _p2 | particle 2 |
[in] | _ecca | eccentricity anomaly |
[in] | _G | gravitational constant |
◆ calcSemiFromPeriod()
void COMM::Binary::calcSemiFromPeriod |
( |
const Float & |
_G | ) |
|
|
inline |
from period calculate semi-major axis
- Parameters
-
[in] | _G | gravitational constant |
◆ evolve()
void COMM::Binary::evolve |
( |
const Float |
_dt | ) |
|
|
inline |
Solve kepler motion for dt.
Evolve current orbit for dt
- Parameters
-
◆ orbitToParticle()
template<class Tptcl >
static void COMM::Binary::orbitToParticle |
( |
Tptcl & |
_p1, |
|
|
Tptcl & |
_p2, |
|
|
const Binary & |
_bin, |
|
|
const Float & |
_ecca, |
|
|
const Float |
_G |
|
) |
| |
|
inlinestatic |
Orbit to position and velocity.
refer to the P3T code developed by Iwasawa M.
- Parameters
-
[out] | _p1 | particle 1 |
[out] | _p2 | particle 2 |
[in] | _bin | binary parameter |
[in] | _ecca | eccentric anomaly (-pi, pi) |
[in] | _G | gravitational constant |
◆ particleToOrbit()
template<class Tptcl >
static void COMM::Binary::particleToOrbit |
( |
Binary & |
_bin, |
|
|
const Tptcl & |
_p1, |
|
|
const Tptcl & |
_p2, |
|
|
const Float |
_G |
|
) |
| |
|
inlinestatic |
position velocity to orbit
◆ particleToSemiEcc()
template<class Tpi , class Tpj >
static void COMM::Binary::particleToSemiEcc |
( |
Float & |
_semi, |
|
|
Float & |
_ecc, |
|
|
Float & |
_r, |
|
|
Float & |
_rv, |
|
|
const Tpi & |
_p1, |
|
|
const Tpj & |
_p2, |
|
|
const Float |
_G |
|
) |
| |
|
inlinestatic |
position velocity to orbit semi-major axis and eccentricity
- Parameters
-
[out] | _semi | semi-major axis |
[out] | _ecc | eccentricity |
[out] | _r | distance between two particles |
[out] | _rv | relative position dot velocity |
[in] | _p1 | particle 1 |
[in] | _p2 | particle 2 |
[in] | _G | gravitational constant |
◆ particleToSemiEccPeriod()
template<class Tpi , class Tpj >
void COMM::Binary::particleToSemiEccPeriod |
( |
const Tpi & |
_p1, |
|
|
const Tpj & |
_p2, |
|
|
const Float |
_G |
|
) |
| |
|
inline |
calcualte semi, ecc and period
- Parameters
-
[in,out] | _p1 | particle 1 |
[in,out] | _p2 | particle 2 |
[in] | _G | gravitational constant |
◆ periodToSemi()
static Float COMM::Binary::periodToSemi |
( |
const Float & |
_period, |
|
|
const Float & |
_mtot, |
|
|
const Float & |
_G |
|
) |
| |
|
inlinestatic |
from period calculate semi-major axis
- Parameters
-
[in] | _period | period |
[in] | _mtot | total mass of binary |
[in] | _G | gravitational constant |
◆ print()
void COMM::Binary::print |
( |
std::ostream & |
_os | ) |
const |
|
inline |
◆ printColumn()
void COMM::Binary::printColumn |
( |
std::ostream & |
_fout, |
|
|
const int |
_width = 20 |
|
) |
| |
|
inline |
print data of class members using column style
print data of class members in one line for column style. Notice no newline is printed at the end
- Parameters
-
[out] | _fout | std::ostream output object |
[in] | _width | print width (defaulted 20) |
◆ printColumnTitle()
static void COMM::Binary::printColumnTitle |
( |
std::ostream & |
_fout, |
|
|
const int |
_width = 20 |
|
) |
| |
|
inlinestatic |
print titles of class members using column style
print titles of class members in one line for column style
- Parameters
-
[out] | _fout | std::ostream output object |
[in] | _width | print width (defaulted 20) |
◆ readAscii()
void COMM::Binary::readAscii |
( |
std::istream & |
_fin | ) |
|
|
inline |
read class data to file with ASCII format
- Parameters
-
[in] | _fin | std::istream file for input |
◆ readBinary()
void COMM::Binary::readBinary |
( |
FILE * |
_fin | ) |
|
|
inline |
read class data to file with binary format
- Parameters
-
[in] | _fin | FILE type file for reading |
◆ rotateToOriginalFrame()
void COMM::Binary::rotateToOriginalFrame |
( |
Float * |
_vec | ) |
|
|
inline |
rotate position vector from binary rest-frame to original frame based on three angles
- Parameters
-
[out] | _vec | vector to rotate |
◆ semiToPeriod()
static Float COMM::Binary::semiToPeriod |
( |
const Float & |
_semi, |
|
|
const Float & |
_mtot, |
|
|
const Float & |
_G |
|
) |
| |
|
inlinestatic |
from semi-major axis calculate period;
- Parameters
-
[in] | _semi | semi-major axis |
[in] | _mtot | total mass of binary |
[in] | _G | gravitational constant |
◆ solveKepler()
static void COMM::Binary::solveKepler |
( |
Binary & |
_bin, |
|
|
const Float |
_dt |
|
) |
| |
|
inlinestatic |
solve kepler orbit after dt
refer to the P3T code developed by Iwasawa M.
- Parameters
-
[in,out] | _bin | binary orbit |
[in] | _dt | evolution time |
◆ writeAscii()
void COMM::Binary::writeAscii |
( |
std::ostream & |
_fout | ) |
const |
|
inline |
write class data to file with ASCII format
- Parameters
-
[in] | _fout | std:osteram file for output |
◆ writeBinary()
void COMM::Binary::writeBinary |
( |
FILE * |
_fp | ) |
const |
|
inline |
write class data to file with binary format
- Parameters
-
[in] | _fp | FILE type file for output |
◆ am
◆ ecc
◆ ecca
◆ incline
Float COMM::Binary::incline |
◆ m1
◆ m2
◆ period
Float COMM::Binary::period |
◆ rot_horizon
Float COMM::Binary::rot_horizon |
◆ rot_self
Float COMM::Binary::rot_self |
◆ semi
◆ stab
◆ t_peri
Float COMM::Binary::t_peri |
The documentation for this class was generated from the following file: