Describes a system of two bodies orbiting each other. More...
#include <BinarySystem.h>
Public Member Functions | |
BinarySystem (DissipatingBody &body1, DissipatingBody &body2, const std::string &system_name="") | |
Construct a binary system. More... | |
const std::string | get_name () const |
Returns the name of the system. More... | |
virtual int | configure (bool initialize, double age, double semimajor, double eccentricity, const double *spin_angmom, const double *inclination, const double *periapsis, Core::EvolModeType evolution_mode) |
Sets the current state of the system. More... | |
int | configure (bool initialize, double age, const double *parameters, Core::EvolModeType evolution_mode) |
double | age () const |
Returns the present age of the system in Gyr. More... | |
const DissipatingBody & | primary () const |
Returns the primary body in the system (const). More... | |
const DissipatingBody & | secondary () const |
Returns the secondary body in the system (const). More... | |
unsigned | number_zones () const |
The total number of zones in both system bodies. More... | |
unsigned | number_locked_zones () const |
How many zones on either body are currently locked. More... | |
double | semimajor () const |
The current semimajor axis of the system. More... | |
double | eccentricity () const |
The current eccentricity of the system. More... | |
Core::EvolModeType | fill_orbit (std::valarray< double > &orbit) const |
Fills an array with the parameters expected by differential_equations() and jacobian(), returning the evolution mode. More... | |
double | above_lock_fraction (unsigned locked_zone_index, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, unsigned deriv_zone_index=0, bool secondary_radius=false) |
The fraction of an infinitesimal timestep that a zone spends spinning faster than the lock it is in. More... | |
int | differential_equations (double age, const double *parameters, Core::EvolModeType evolution_mode, double *differential_equations, bool expansion_error=false) |
The differential equation and jacobian for the evolution of the system. More... | |
int | jacobian (double age, const double *parameters, Core::EvolModeType evolution_mode, double *param_derivs, double *age_derivs) |
void | check_for_lock (int orbital_freq_mult, int spin_freq_mult, unsigned short body_index, unsigned zone_index, short direction) |
Check if a spin-orbit lock can be held and updates the system as necessary to calculate subsequent evolution. More... | |
virtual double | minimum_separation (bool deriv=false) const |
Smallest semimajor axis at which the secondary can survive for the latest system configuration. More... | |
Core::EvolModeType | evolution_mode () |
The evolution mode of last call to configure(). More... | |
virtual void | secondary_died () |
Update the system to account for the death of the secondary. More... | |
virtual void | release_lock (unsigned locked_zone_index, short direction) |
Releases the lock to one of the locked zones. More... | |
virtual void | add_to_evolution () |
Appends the state defined by last configure(), to the evolution. More... | |
virtual void | reset_evolution () |
Resets the evolution of the system. More... | |
virtual void | rewind_evolution (unsigned nsteps) |
Discards the last steps from the evolution. More... | |
virtual CombinedStoppingCondition * | stopping_conditions () |
Conditions detecting the next possible doscontinuity in the evolution. More... | |
virtual void | reached_critical_age (double age) |
Change the system as necessary at the given age. More... | |
virtual double | next_stop_age () const |
The next age when the evolution needs to be stopped for a system change. More... | |
const std::list< double > & | semimajor_evolution () const |
The tabulated evolution of the semimajor axis so far. More... | |
const std::list< double > & | semimajor_evolution_rate () const |
The tabulated evolution of the semimajor axis so far. More... | |
const std::list< double > & | eccentricity_evolution () const |
The tabulated evolution of the eccentricity so far. More... | |
const std::list< double > & | eccentricity_evolution_rate () const |
The tabulated evolution of the eccentricity so far. More... | |
virtual void | change_e_order (unsigned new_e_order) |
Change the eccentricity expansion order for all dissipative zones. More... | |
virtual unsigned | eccentricity_order () const |
Private Member Functions | |
void | find_locked_zones () |
Fills the __locked_zones list. More... | |
int | locked_surface_differential_equations (double *evolution_rates, bool expansion_error) const |
Differential equations for the rotation of the zones of body 0 with the topmost zone rotating with a fixed frequency. More... | |
void | locked_surface_jacobian (double *param_derivs, double *age_derivs) const |
Jacobian for the evolution of the rotation of the zones of body 0 with the topmost zone rotating with a fixed frequency. More... | |
int | single_body_differential_equations (double *evolution_rates, bool expansion_error) const |
Differential equations for the rotation of the zones of body 0 if no other body is present. More... | |
void | fill_single_body_jacobian (double *inclination_param_derivs, double *periapsis_param_derivs, double *angmom_param_derivs, double *inclination_age_derivs, double *periapsis_age_derivs, double *angmom_age_derivs) const |
Fills the jacobian for a system consisting of one isolated body. More... | |
void | single_body_jacobian (double *param_derivs, double *age_derivs) const |
Jacobian for the evolution of the rotation of the zones of body 0 if no other body is present. More... | |
double | semimajor_evolution (double orbit_power, double orbit_power_deriv=Core::NaN) const |
Returns the rate of evolution of the semimajor axis or one of its derivatives. More... | |
double | semimajor_evolution_expansion_error () const |
Estimate of the error in the value returned by semimajor_evolution() due to truncating the tidal potential eccentricity expansion. More... | |
double | eccentricity_evolution (double orbit_power, double orbit_angmom_gain, double orbit_power_deriv=Core::NaN, double orbit_angmom_gain_deriv=Core::NaN, bool semimajor_deriv=true) const |
Returns the rate of evolution of the eccentricity or one of its derivatives. More... | |
double | eccentricity_evolution_expansion_error () const |
Estimate of the error in the value returned by eccentricity_evolution() due to truncating the tidal potential eccentricity expansion. More... | |
void | above_lock_problem_deriv_correction (Dissipation::QuantityEntry entry, bool body1_deriv, Eigen::MatrixXd &matrix, Eigen::VectorXd &rhs) const |
Makes corrections to the matrix and RHS to accomodate the given derivative for the linear problem that defines the above fractions. More... | |
void | calculate_above_lock_fractions (Eigen::VectorXd &fractions, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, bool body1_deriv=true) |
Calculates the fraction of a timestep above spin-orbit lock for all locked zones. More... | |
Eigen::VectorXd | above_lock_fractions_deriv (Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_index) |
Calculates derivatives of the above lock fractions w.r.t. quantities of non-surface zones. More... | |
void | fill_above_lock_fractions_deriv () |
Fills the __above_lock_fractinos_*_deriv members. More... | |
void | update_above_lock_fractions () |
Solves for and sets the above lock fractions and their derivatives. More... | |
void | fill_orbit_torque_and_power () |
Update the values of ::__orbit_power, ::__orbit_torque, ::__orbit_angmom_gain and their associated expansion errors. More... | |
int | binary_differential_equations (double *differential_equations, bool expansion_error) const |
The differential equations for a system with both bodies present. More... | |
template<typename VALUE_TYPE > | |
void | add_body_rate_deriv (const DissipatingBody &body, VALUE_TYPE(DissipatingBody::*func)(Dissipation::QuantityEntry, unsigned, const Eigen::VectorXd &) const, std::valarray< VALUE_TYPE > &orbit_rate_deriv, unsigned offset) const |
Adds the derivatives of a rate by which the orbit is changing due to tidal interactions with a single body to an array. More... | |
void | fill_orbit_power_deriv (std::valarray< double > &orbit_power_deriv) const |
Computes the derivatives w.r.t. the evolution quantities of the orbit energy gain. More... | |
void | fill_orbit_angmom_gain_deriv (std::valarray< double > &orbit_angmom_gain_deriv) const |
Computes the derivatives w.r.t. the evolution quantities of the orbit angular momentum gain. More... | |
void | semimajor_jacobian (const std::valarray< double > &orbit_power_deriv, bool a6p5, double *param_derivs, double &age_deriv) const |
Computes the row of the jacobian corresponding to the semimajor axis. More... | |
void | eccentricity_jacobian (const std::valarray< double > &orbit_power_deriv, const std::valarray< double > &orbit_angmom_gain_deriv, bool a6p5, double *param_derivs, double &age_deriv) const |
Computes the row of the jacobian corresponding to the eccentricity. More... | |
void | angle_evolution_age_deriv (DissipatingBody &body, unsigned zone_ind, double sin_inc, double cos_inc, unsigned locked_zone_ind, double &inclination, double &periapsis) const |
Computes the partial age derivatives of the inclination or periapsis evolutiono equations. More... | |
void | angle_evolution_orbit_deriv (Dissipation::QuantityEntry entry, double angmom_deriv, DissipatingBody &body, unsigned zone_ind, double sin_inc, double cos_inc, unsigned locked_zone_ind, double &inclination, double &periapsis) const |
Computes the partial derivatives w.r.t. semimamjor axis or eccentricity of the inclination and periapsis evolution equations. More... | |
void | fill_orbit_torque_deriv (Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_ind, std::valarray< Eigen::Vector3d > &orbit_torque_deriv) const |
void | fill_zone_torque_deriv (Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_ind, std::valarray< Eigen::Vector3d > &zone_torque_deriv) const |
Calculates the derivatives of the torque on a zone w.r.t. a zone-specific quantity. More... | |
void | inclination_evolution_zone_derivs (Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_ind, double zone_x_torque_above, double zone_x_torque_below, const std::valarray< Eigen::Vector3d > &zone_torque_deriv, const Eigen::Vector3d &orbit_torque, const std::valarray< Eigen::Vector3d > &orbit_torque_deriv, const std::valarray< Eigen::VectorXd > &above_frac_deriv, double sin_inc, double cos_inc, unsigned locked_zone_ind, double *result) const |
Calculates the derivatives of an inclination evolution equation w.r.t. a zone specific quentity. More... | |
void | periapsis_evolution_zone_derivs (Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_ind, double zone_y_torque_above, double zone_y_torque_below, const std::valarray< Eigen::Vector3d > &zone_torque_deriv, double orbit_y_torque, const std::valarray< Eigen::Vector3d > &orbit_torque_deriv, const std::valarray< Eigen::VectorXd > &above_frac_deriv, double sin_inc, double cos_inc, unsigned locked_zone_ind, double *result) const |
Computes the derivatives of the periapsis evolution w.r.t. a zone specific quantity for all zones. More... | |
void | spin_angmom_evolution_zone_derivs (Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_ind, double zone_z_torque_above, double zone_z_torque_below, const std::valarray< Eigen::Vector3d > &zone_torque_deriv, const std::valarray< Eigen::VectorXd > &above_frac_deriv, unsigned locked_zone_ind, double *result) const |
void | binary_jacobian (double *param_derivs, double *age_derivs) const |
Calculates the jacobian for the evolution of a system of two bodies. More... | |
void | fill_locked_surface_orbit (std::valarray< double > &orbit) const |
Implements fill_orbit() for LOCKED_SURFACE_SPIN evolution mode. More... | |
void | fill_binary_orbit (std::valarray< double > &orbit) const |
Implements fill_orbit() for BINARY evolution mode. More... | |
void | fill_single_orbit (std::valarray< double > &orbit) const |
Implements fill_orbit() for SINGLE evolution mode. More... | |
Private Attributes | |
std::string | __name |
The name of the binary system (e.g. "HAT-P-20") More... | |
std::list< double > | __semimajor_evolution |
The evolution of the semimajor axis recorded by add_to_evolution() so far. More... | |
std::list< double > | __eccentricity_evolution |
The evolution of the eccentricity recorded by add_to_evolution() so far. More... | |
std::list< double > | __semimajor_rate_evolution |
The rate at which the semimajor axis evolved at each. More... | |
std::list< double > | __eccentricity_rate_evolution |
The rate at which the eccentricity evolved at each. More... | |
double | __age |
The present age of the stellar system in Gyrs. More... | |
double | __semimajor |
The current semimajor axis. More... | |
double | __eccentricity |
The current eccentricity. More... | |
double | __orbital_energy |
The current orbital energy. More... | |
double | __orbital_angmom |
The current orbital angular momentum. More... | |
double | __orbit_power |
The rate at which the orbit gains energy due to tides. More... | |
double | __orbit_power_expansion_error |
Estimate of the error in __orbit_power due to truncating the tidal potential eccentricity expansion. More... | |
double | __orbit_angmom_gain |
The rate at which the orbit gains angular momentum due to tides. More... | |
double | __orbit_angmom_gain_expansion_error |
Estimate of the error in __orbit_angmom_gain due to truncating the tidal potential eccentricity expansion. More... | |
double | __semimajor_rate |
The current rate at which the semiamjor axis is changing. More... | |
double | __eccentricity_rate |
The current rate at which the eccentricity is changing. More... | |
Eigen::Vector3d | __orbit_torque |
The torque on the orbit in the coordinate system of the outermost zone of the first body. More... | |
Eigen::Vector3d | __orbit_torque_expansion_error |
An estiamte of the error in ::__orbit_torque due to truncating the eccentricity series of the tidal potential. More... | |
Core::EvolModeType | __evolution_mode |
The evolution mode from the last call to configure();. More... | |
std::list< unsigned > | __locked_zones |
A list of indices of locked zones. More... | |
std::valarray< Eigen::VectorXd > | __above_lock_fractions |
The above lock fractinos for the locked zones and their derivatives. More... | |
std::valarray< Eigen::VectorXd > | __above_lock_fractions_inclination_deriv |
The derivatives of the above lock fractions w.r.t. the inclinations of the zones. More... | |
std::valarray< Eigen::VectorXd > | __above_lock_fractions_periapsis_deriv |
The derivatives of the above lock fractions w.r.t. the periapses of the zones. More... | |
std::valarray< Eigen::VectorXd > | __above_lock_fractions_inertia_deriv |
The derivatives of the above lock fractions w.r.t. the moments of inertia of the zones. More... | |
std::valarray< Eigen::VectorXd > | __above_lock_fractions_angmom_deriv |
The derivatives of the above lock fractions w.r.t. the angular momenta of the zones. More... | |
Eigen::VectorXd | __above_lock_fractions_body2_radius_deriv |
The derivative of the above lock fractions w.r.t. to the radius of the secondardy. More... | |
Eigen::ColPivHouseholderQR< Eigen::MatrixXd > | __above_lock_fractions_decomp |
The matrix that defines the problem to solve for __above_lock_fractions. More... | |
DissipatingBody & | __body1 |
The first body in the system. More... | |
DissipatingBody & | __body2 |
The second body in the system. More... | |
Describes a system of two bodies orbiting each other.
The following variable are referred to throughout:
The index of the outermost zone of a body is zero and increases inward.
All rates of change are per Gyr
Positive tidal power and tidal torque mean that energy is being removed from the orbit and deposited into the bodies.
Definition at line 56 of file BinarySystem.h.
|
inline |
Construct a binary system.
body1 | The first body in the system. Assumed to always be there, so for a star-planet system this should be the star. |
body2 | The second body in the system, initially may not be there and later may be engulfed by the first body. |
system_name | The name of the system. |
Definition at line 763 of file BinarySystem.h.
double Evolve::BinarySystem::above_lock_fraction | ( | unsigned | locked_zone_index, |
Dissipation::QuantityEntry | entry = Dissipation::NO_DERIV , |
||
unsigned | deriv_zone_index = 0 , |
||
bool | secondary_radius = false |
||
) |
The fraction of an infinitesimal timestep that a zone spends spinning faster than the lock it is in.
locked_zone_index | The index within the list of locked zones of the zone we need. |
entry | What derivative to return. Only non-zone specific quantities and, Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::SPIN_ANGMOM and Dissipation::MOMENT_OF_INERTIA are allowed. |
deriv_zone_index | The index of the zone whose quantity we want the derivative w.r.t. for zone specific quantities. |
secondary_radius | If deriv==Dissipation::RADIUS this argument determines the body whose radius the derivative is w.r.t. |
Definition at line 2073 of file BinarySystem.cpp.
|
private |
Calculates derivatives of the above lock fractions w.r.t. quantities of non-surface zones.
entry | The derivative to calculate. Only INCLINATION, PERIAPSIS, MOMENT_OF_INERTIA and SPIN_ANGMOM derivatives are supported. |
body | The body whose zone's inclination we want the derivative w.r.t. |
zone_index | The zone whose inclination we are finding the derivative w.r.t. |
Definition at line 700 of file BinarySystem.cpp.
|
private |
Makes corrections to the matrix and RHS to accomodate the given derivative for the linear problem that defines the above fractions.
entry | The derivative being calculated |
body1_deriv | For zone-specific quantities, are we differentiating w.r.t. a zone part of body1? |
matrix | The matrix to update. |
rhs | The RHS vector to update. |
Definition at line 501 of file BinarySystem.cpp.
|
private |
Adds the derivatives of a rate by which the orbit is changing due to tidal interactions with a single body to an array.
With the appropriate arguments this can either add the rate of change of the orbital energy or the torque on the orbit due to the tides on one body to the given array.
body | The body whose contribution we wish to add. |
func | The quantity we are collecting. Should be either body.tidal_orbit_power or body.tidal_orbit_torque. |
orbit_rate_deriv | The array to add to. The indices match the parameters being evolved. |
offset | If body is __body1 this should be zero, if it is __body2, it should be the number of zones in __body1. |
Definition at line 1047 of file BinarySystem.cpp.
|
virtual |
Appends the state defined by last configure(), to the evolution.
Definition at line 2353 of file BinarySystem.cpp.
|
inline |
Returns the present age of the system in Gyr.
Definition at line 833 of file BinarySystem.h.
|
private |
Computes the partial age derivatives of the inclination or periapsis evolutiono equations.
body | The body whose zone to calculate the age derivative for. |
zone_ind | The index of the zone in body to calculate the age derivative for. |
sin_inc | The sin of the inclination between the zone and the orbit. |
cos_inc | The cos of the inclination between the zone and the orbit. |
locked_zone_ind | The index of the current zone in the list of locked zones (used only if zone is locked, ignored otherwise) |
inclination | Destination overwritten by the age derivative of the inclination evolution rate. |
periapsis | Destination overwritten by the age derivative of the periapsis evolution rate. |
Definition at line 1209 of file BinarySystem.cpp.
|
private |
Computes the partial derivatives w.r.t. semimamjor axis or eccentricity of the inclination and periapsis evolution equations.
entry | The derivative to compute, should be either Dissipation::SEMIMAJOR or Dissipation::ECCENTRICITY |
angmom_deriv | The derivative of the orbital angular momentum w.r.t. deriv. |
body | The body whose zone's evolution we are differentiating. |
zone_ind | The index of the zone whose evolution we are differentiating. |
sin_inc | The sin of the inclination between the zone and the orbit. |
cos_inc | The cos of the inclination between the zone and the orbit. |
locked_zone_ind | The index of the current zone in the list of locked zones (used only if zone is locked, ignored otherwise) |
inclination | Overwritten by the derivative of the inclination evolution rate. |
periapsis | Overwritten by the derivative of the periapsis evolution rate. |
Definition at line 1287 of file BinarySystem.cpp.
|
private |
The differential equations for a system with both bodies present.
Also updates the current evolution rate private members if expansion_error is false.
differential_equations | On output is set to the rates of change of the evolution variables. See differintal_equations() for details. |
expansion_error | If true, instead of returning the evolution rates, returns an estimete of the error in those due to truncating the eccentricity expansion series of the tidal potential. |
Definition at line 895 of file BinarySystem.cpp.
|
private |
Calculates the jacobian for the evolution of a system of two bodies.
param_derivs | On output is set to the Jacobian. |
age_derivs | On output is set to the partial age derivatives of the evolution equations. |
Definition at line 1679 of file BinarySystem.cpp.
|
private |
Calculates the fraction of a timestep above spin-orbit lock for all locked zones.
The configure() method must already have been called.
Only locked zones get an entry. The locked zones of body1 are first from outside to inside, followed by the locked zones for body2.
fractions | The vector to fill with the calculated fractions. |
entry | The derivative of the above lock fractions to calculate. For zone-dependent quantities the derivatives are w.r.t. the surface zone quantity. |
body1_deriv | Only used if deriv indicates a zone-dependent quantity. In that case, it specifies the body whose surface zone we are differentiating against. |
Definition at line 597 of file BinarySystem.cpp.
|
inlinevirtual |
Change the eccentricity expansion order for all dissipative zones.
new_e_order | The new eccentricity expansion order. |
Definition at line 1060 of file BinarySystem.h.
void Evolve::BinarySystem::check_for_lock | ( | int | orbital_freq_mult, |
int | spin_freq_mult, | ||
unsigned short | body_index, | ||
unsigned | zone_index, | ||
short | direction | ||
) |
Check if a spin-orbit lock can be held and updates the system as necessary to calculate subsequent evolution.
orbital_freq_mult | The multiplier of the orbital frequency at the lock |
spin_freq_mult | The multiplier of the spin frequency of the potentially locked zone at the lock. |
body_index | The body whose zone is being locked. |
zone_index | The zone being unlocked. |
direction | The direction in which the lock will be left if it does not hold. |
Definition at line 2169 of file BinarySystem.cpp.
|
virtual |
Sets the current state of the system.
initialize | Is this the first time configure() is invoked? |
age | The age to set the system to. |
semimajor | The semimajor axis of the orbit. |
eccentricity | The eccentricity of the orbit. |
spin_angmom | The spin angular momenta of the zones of the bodies (body 1 first, outermost zone to innermost, followed by body 2). |
inclination | The inclinations of the zones of the bodies (same order as spin_angmom). The surface zone inclination must be omitted for single body systems. |
periapsis | The arguments of periapsis of the zones of the bodies (same order as spin_angmom, but not including the surface zone of the first body). |
evolution_mode | The evolution mode to assume. |
Definition at line 1912 of file BinarySystem.cpp.
int Evolve::BinarySystem::configure | ( | bool | initialize, |
double | age, | ||
const double * | parameters, | ||
Core::EvolModeType | evolution_mode | ||
) |
Sets the current state of the system directly from the evolutino variables.
initialize | Is this the first time configure() is invoked? |
age | The current age of the system. |
parameters | The parameters being evolved. See differential_equations(). |
evolution_mode | The evolution mode to assume for the system. |
Definition at line 2005 of file BinarySystem.cpp.
int Evolve::BinarySystem::differential_equations | ( | double | age, |
const double * | parameters, | ||
Core::EvolModeType | evolution_mode, | ||
double * | differential_equations, | ||
bool | expansion_error = false |
||
) |
The differential equation and jacobian for the evolution of the system.
Deliberate side effects:
age | The current age of the system. |
parameters | Contains the variables being evolved. The expected content depends on the values of evolution_mode and mhether any spin-orbit locks are held. |
If evolution_mode is LOCKED_SURFACE_SPIN:
If evolution mode is SINGLE
The TABULATION evolution mode is not allowed.
If evolution mode is BINARY
evolution_mode | The evolution mode to assume for the system. |
differential_equations | On outputs gets filled with the rates at which the entries in parameters evolve. It is assumed that sufficient space has been allocated to hold the results. |
expansion_error | If true, instead of returning the evolution rates, returns an estimete of the error in those due to truncating the eccentricity expansion series of the tidal potential. |
Definition at line 2106 of file BinarySystem.cpp.
|
inline |
The current eccentricity of the system.
Definition at line 855 of file BinarySystem.h.
|
private |
Returns the rate of evolution of the eccentricity or one of its derivatives.
The orbit and age of the system must already be set by calling configure().
Supports differentiation with repesct to the semimajor axis and eccentricity. If derivatives with respect to other quantities are required, simply provide the derivative of the tidal power and torque with respect to the quantity as the first and second arguments.
orbit_power | See semimajor_evolution() |
orbit_angmom_gain | The rate at which the orbit gains angular momentum (total for all zones of all bodies) in \(M_\odot R_\odot^2 \mathrm{day}^{-2}\mathrm{Gyr}^{-1}\) |
orbit_power_deriv | If this is not NaN, the derivative of the eccentricity evolution rate w.r.t. either the semimajor axis or the case, this value must be the derivative of orbit_power w.r.t. the same this as the desired derivative. |
orbit_angmom_gain_deriv | If orbit_power_deriv is not NaN, this must be set to the derivative of orbit_torque with respect to the same variable as orbit_power_deriv. |
semimajor_deriv | If true the derivative calculated is assumed to be w.r.t. the semimajor axis. |
Definition at line 430 of file BinarySystem.cpp.
|
inline |
The tabulated evolution of the eccentricity so far.
Definition at line 1052 of file BinarySystem.h.
|
private |
Estimate of the error in the value returned by eccentricity_evolution() due to truncating the tidal potential eccentricity expansion.
Definition at line 485 of file BinarySystem.cpp.
|
inline |
The tabulated evolution of the eccentricity so far.
Definition at line 1056 of file BinarySystem.h.
|
private |
Computes the row of the jacobian corresponding to the eccentricity.
orbit_power_deriv | The derivatives of the orbit energy gain w.r.t. the evolution variables and age (last entry). |
orbit_angmom_gain_deriv | The derivatives of the orbit angular momentum gain w.r.t. the evolution variables and age (last entry). |
a6p5 | Is the first variable \(a^{6.5}\) instead of a? |
param_derivs | The location to fill with the deriavtives of the eccentricity evolution equation w.r.t. the evolution variables. |
age_deriv | The location to set to the age derivative of the eccentricity evolution equation. |
Definition at line 1182 of file BinarySystem.cpp.
|
virtual |
To what order should eccentricity expansion be performed for the given value of the eccentricity.
Definition at line 2408 of file BinarySystem.cpp.
|
inline |
The evolution mode of last call to configure().
Definition at line 996 of file BinarySystem.h.
|
private |
Fills the __above_lock_fractinos_*_deriv members.
Definition at line 779 of file BinarySystem.cpp.
|
private |
Implements fill_orbit() for BINARY evolution mode.
Definition at line 1855 of file BinarySystem.cpp.
|
private |
Implements fill_orbit() for LOCKED_SURFACE_SPIN evolution mode.
Definition at line 1839 of file BinarySystem.cpp.
Core::EvolModeType Evolve::BinarySystem::fill_orbit | ( | std::valarray< double > & | orbit | ) | const |
Fills an array with the parameters expected by differential_equations() and jacobian(), returning the evolution mode.
The system must be appropriately configure() -ed already.
orbit | The orbit to fill (resized as necessary). |
Definition at line 2062 of file BinarySystem.cpp.
|
private |
Computes the derivatives w.r.t. the evolution quantities of the orbit angular momentum gain.
orbit_angmom_gain_deriv | Location to fill. |
Definition at line 1122 of file BinarySystem.cpp.
|
private |
Computes the derivatives w.r.t. the evolution quantities of the orbit energy gain.
orbit_power_deriv | Location to fill. |
Definition at line 1107 of file BinarySystem.cpp.
|
private |
Update the values of ::__orbit_power, ::__orbit_torque, ::__orbit_angmom_gain and their associated expansion errors.
Called as part of ::configure().
Definition at line 849 of file BinarySystem.cpp.
|
private |
entry | The derivatives to compute. Sholud be one of: Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::MOMENT_OF_INERTIA, Dissipation::SPIN_ANGMOM |
body | The body whose zone's coordinate system we are expressing the torque in. |
zone_ind | The index of the zone whose coordinate system we are expressing the torque in. |
orbit_torque_deriv | On output gets filled with the derivatives of the torque on the orbit in the coordinate system of the zone identified by body and zone_ind w.r.t. deriv. Should be indexed by zone, with __body1 zones first. |
Definition at line 1357 of file BinarySystem.cpp.
|
private |
Fills the jacobian for a system consisting of one isolated body.
inclination_param_derivs | The rows of the jacobian corresponding to the zone inclinations. |
periapsis_param_derivs | The rows of the jacobian corresponding to the zone periapses. |
angmom_param_derivs | The rows of the jacobian corresponding to the zone angular momenta. |
inclination_age_derivs | The part of the age derivatives array corresponding to the zone inclinations. |
periapsis_age_derivs | The part of the age derivatives array corresponding to the zone periapses. |
angmom_age_derivs | The part of the age derivatives array corresponding to the zone angular momenta. |
Definition at line 229 of file BinarySystem.cpp.
|
private |
Implements fill_orbit() for SINGLE evolution mode.
Definition at line 1886 of file BinarySystem.cpp.
|
private |
Calculates the derivatives of the torque on a zone w.r.t. a zone-specific quantity.
entry | The derivatives to compute. Sholud be one of: Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::MOMENT_OF_INERTIA, Dissipation::SPIN_ANGMOM |
body | The body whose zone's evolution we are differentiating. |
zone_ind | The index of the zone whose evolution we are differentiating. |
zone_torque_deriv | Overwritten with the result. If the zone is not locked this contains the derivatives w.r.t. to the quantity for the zone above, the zone and the zone below. If the zone is locked, the first three quantities are for the torque below the lock and an additional quantity is added at the end givin the derivative of the torque above the lock w.r.t. to the quantity of this zone. |
Definition at line 1430 of file BinarySystem.cpp.
|
private |
Fills the __locked_zones list.
Definition at line 13 of file BinarySystem.cpp.
|
inline |
Returns the name of the system.
Definition at line 782 of file BinarySystem.h.
|
private |
Calculates the derivatives of an inclination evolution equation w.r.t. a zone specific quentity.
entry | The derivatives to compute. Sholud be one of: Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::MOMENT_OF_INERTIA, Dissipation::SPIN_ANGMOM |
body | The body whose zone's evolution we are differentiating. |
zone_ind | The index of the zone whose evolution we are differentiating. |
zone_x_torque_above | The x component of the torque on the zone (above a potential spin orbit lock). Ignored if the zone is not locked. |
zone_x_torque_below | The x component of the torque on the zone (below a potential spin orbit lock). |
zone_torque_deriv | To be computed by fill_zone_torque_deriv(). |
orbit_torque | The torque on the orbit in the coordinate system of the zone defined by the above 2 arguments. |
orbit_torque_deriv | The result of fill_orbit_torque_deriv. |
above_frac_deriv | The derivatives of the above_lock fractions. |
sin_inc | The sin of the inclination between the zone and the orbit. |
cos_inc | The cos of the inclination between the zone and the orbit. |
locked_zone_ind | The index of the current zone in the list of locked zones (used only if zone is locked, ignored otherwise) |
result | Overwritten by the result. |
Definition at line 1462 of file BinarySystem.cpp.
int Evolve::BinarySystem::jacobian | ( | double | age, |
const double * | parameters, | ||
Core::EvolModeType | evolution_mode, | ||
double * | param_derivs, | ||
double * | age_derivs | ||
) |
The jacobian of the evolution equations.
Calls configure()!
age | The age of the system. |
parameters | See differential_equations() |
evolution_mode | The evolution mode to assume for the system. |
param_derivs | The matrix of partial derivatives of the evolution rates for the parameters w.r.t. each of the parameters. |
age_derivs | The partialderivatives of the evolution rates for the parameters w.r.t. age. |
Definition at line 2146 of file BinarySystem.cpp.
|
private |
Differential equations for the rotation of the zones of body 0 with the topmost zone rotating with a fixed frequency.
The age of the system must already be set appropriately by configure().
evolution_rates | On output is set to the rates of change of \(S^0_i\). |
expansion_error | If true, instead of returning the evolution rates, returns an estimete of the error in those due to truncating the eccentricity expansion series of the tidal potential. |
Definition at line 31 of file BinarySystem.cpp.
|
private |
Jacobian for the evolution of the rotation of the zones of body 0 with the topmost zone rotating with a fixed frequency.
The age of the system must already be set appropriately by configure().
param_derivs | On output is set to the Jacobian. |
age_derivs | On output is set to the partial age derivatives of the evolution equations. |
Definition at line 68 of file BinarySystem.cpp.
|
virtual |
Smallest semimajor axis at which the secondary can survive for the latest system configuration.
By default returns the larger of:
deriv | If true the rate of change (per Gyr) is returned. |
Definition at line 2248 of file BinarySystem.cpp.
|
virtual |
The next age when the evolution needs to be stopped for a system change.
Reimplemented in Evolve::DiskBinarySystem.
Definition at line 2402 of file BinarySystem.cpp.
|
inline |
How many zones on either body are currently locked.
Definition at line 848 of file BinarySystem.h.
|
inline |
The total number of zones in both system bodies.
Definition at line 842 of file BinarySystem.h.
|
private |
Computes the derivatives of the periapsis evolution w.r.t. a zone specific quantity for all zones.
entry | The derivative to compute. |
body | The body whose zone's periapsis evolution to differentiate. |
zone_ind | The zone whose periapsis evolution to differentiate. |
zone_y_torque_above | The y component of the torque on the zone above a lock if locked. |
zone_y_torque_below | The y component of the torque on the zone (assuming below a lock if locked). |
zone_torque_deriv | The derivative of the zone torque. |
orbit_y_torque | The y torque on the orbit. |
orbit_torque_deriv | The derivative of the orbit torque. |
above_frac_deriv | The derivatives of the above lock fractions |
sin_inc | The sin of the inclination between the zone and the orbit. |
cos_inc | The cos of the inclination between the zone and the orbit. |
locked_zone_ind | The index in the list of locked zones of this zone if it is locked. |
result | Overwritten by the result. Should have the correct size. |
Definition at line 1542 of file BinarySystem.cpp.
|
inline |
Returns the primary body in the system (const).
Definition at line 836 of file BinarySystem.h.
|
virtual |
Change the system as necessary at the given age.
Handles things like the disk dissipating, the planet forming and interpolation discontinuities.
Reimplemented in Evolve::DiskBinarySystem.
Definition at line 2396 of file BinarySystem.cpp.
|
virtual |
Releases the lock to one of the locked zones.
locked_zone_index | The index within the list of locked zones of the zone to unlock. |
direction | The direction in which the zone's spin will evolve in the future relative to the lock. |
Definition at line 2321 of file BinarySystem.cpp.
|
virtual |
Resets the evolution of the system.
Definition at line 2363 of file BinarySystem.cpp.
|
virtual |
Discards the last steps from the evolution.
nsteps | How many steps of evolution to discard. |
Definition at line 2373 of file BinarySystem.cpp.
|
inline |
Returns the secondary body in the system (const).
Definition at line 839 of file BinarySystem.h.
|
virtual |
Update the system to account for the death of the secondary.
Adds the entire orbital angular momentum to the surface zone of the primary and enters single body evolution mode.
Definition at line 2271 of file BinarySystem.cpp.
|
inline |
The current semimajor axis of the system.
Definition at line 852 of file BinarySystem.h.
|
private |
Returns the rate of evolution of the semimajor axis or one of its derivatives.
The age and the orbit of the system must already be set by calling configure().
Supports differentiation with repesct to the semimajor axis. If derivatives with respect to other quantities are required, simply provide the derivative of the tidal power with respect to the quantity as the first argument.
orbit_power | The rate at which the orbit gains energy (total for all zones of all bodies) in \(M_\odot R_\odot^2 \mathrm{day}^{-2}\mathrm{Gyr}^{-1}\) |
orbit_power_deriv | If not NaN, the derivative with respect to the semimajoir axis is returned, assuming that this is the derivative of orbit_power with respect to the semimajor axis. |
Definition at line 407 of file BinarySystem.cpp.
|
inline |
The tabulated evolution of the semimajor axis so far.
Definition at line 1044 of file BinarySystem.h.
|
inlineprivate |
Estimate of the error in the value returned by semimajor_evolution() due to truncating the tidal potential eccentricity expansion.
Definition at line 294 of file BinarySystem.h.
|
inline |
The tabulated evolution of the semimajor axis so far.
Definition at line 1048 of file BinarySystem.h.
|
private |
Computes the row of the jacobian corresponding to the semimajor axis.
orbit_power_deriv | The derivatives of the orbit energy gain w.r.t. the evolution variables and age (last entry). |
a6p5 | Is the first variable \(a^{6.5}\) instead of a? |
param_derivs | The location to fill with the deriavtives of the semimajor evolution equation w.r.t. the evolution variables. |
age_deriv | The location to set to the age derivative of the semimajor evolution equation. |
Definition at line 1165 of file BinarySystem.cpp.
|
private |
Differential equations for the rotation of the zones of body 0 if no other body is present.
configure() must already have been called.
evolution_rates | On outputs is set to the rate of change of the orbital parameters. |
expansion_error | If true, instead of returning the evolution rates, returns an estimete of the error in those due to truncating the eccentricity expansion series of the tidal potential. |
Definition at line 163 of file BinarySystem.cpp.
|
private |
Jacobian for the evolution of the rotation of the zones of body 0 if no other body is present.
The age of the system must already be set appropriately by configure().
param_derivs | On output is set to the Jacobian. |
age_derivs | On output is set to the partial age derivatives of the evolution equations. |
Definition at line 394 of file BinarySystem.cpp.
|
private |
entry | The derivative to compute. |
body | The body whose zone's periapsis evolution to differentiate. |
zone_ind | The zone whose periapsis evolution to differentiate. |
zone_z_torque_above | The z component of the torque on the zone above a lock if locked. |
zone_z_torque_below | The z component of the torque on the zone (assuming below a lock if locked). |
zone_torque_deriv | The derivative of the zone torque. |
above_frac_deriv | The derivatives of the above lock fractions |
locked_zone_ind | The index in the list of locked zones of this zone if it is locked. |
result | Overwritten by the result. Should have the correct size. |
Definition at line 1635 of file BinarySystem.cpp.
|
virtual |
Conditions detecting the next possible doscontinuity in the evolution.
Must be deleted when no longer necessary.
Definition at line 2385 of file BinarySystem.cpp.
|
private |
Solves for and sets the above lock fractions and their derivatives.
The system must already be configure() -ed in in BINARY evolution mode.
Definition at line 821 of file BinarySystem.cpp.
|
private |
The above lock fractinos for the locked zones and their derivatives.
For zone-dependent quantities, the derivatives are w.r.t. to the surface zone of __body1.
Definition at line 136 of file BinarySystem.h.
|
private |
The derivatives of the above lock fractions w.r.t. the angular momenta of the zones.
Definition at line 136 of file BinarySystem.h.
|
private |
The derivative of the above lock fractions w.r.t. to the radius of the secondardy.
The derivative w.r.t. to the radius of the primary is already available as __above_lock_fractions[Dissipation::RADIUS].
Definition at line 159 of file BinarySystem.h.
|
private |
The matrix that defines the problem to solve for __above_lock_fractions.
Useful when calculating derivatives.
Definition at line 166 of file BinarySystem.h.
|
private |
The derivatives of the above lock fractions w.r.t. the inclinations of the zones.
Definition at line 136 of file BinarySystem.h.
|
private |
The derivatives of the above lock fractions w.r.t. the moments of inertia of the zones.
Definition at line 136 of file BinarySystem.h.
|
private |
The derivatives of the above lock fractions w.r.t. the periapses of the zones.
Definition at line 136 of file BinarySystem.h.
|
private |
The present age of the stellar system in Gyrs.
Definition at line 80 of file BinarySystem.h.
|
private |
The first body in the system.
This is the body that is allowed to exist both before and after the second body does. It is also the body that should support surface locks.
Definition at line 174 of file BinarySystem.h.
|
private |
The second body in the system.
This body is never allowed to exist alone in the system or to have its surface rotation locked.
Definition at line 174 of file BinarySystem.h.
|
private |
The current eccentricity.
Definition at line 80 of file BinarySystem.h.
|
private |
The evolution of the eccentricity recorded by add_to_evolution() so far.
Definition at line 64 of file BinarySystem.h.
|
mutableprivate |
The current rate at which the eccentricity is changing.
Definition at line 111 of file BinarySystem.h.
|
private |
The rate at which the eccentricity evolved at each.
Definition at line 64 of file BinarySystem.h.
|
private |
The evolution mode from the last call to configure();.
Definition at line 126 of file BinarySystem.h.
|
private |
A list of indices of locked zones.
Definition at line 129 of file BinarySystem.h.
|
private |
The name of the binary system (e.g. "HAT-P-20")
Definition at line 59 of file BinarySystem.h.
|
private |
The rate at which the orbit gains angular momentum due to tides.
Definition at line 80 of file BinarySystem.h.
|
private |
Estimate of the error in __orbit_angmom_gain due to truncating the tidal potential eccentricity expansion.
Definition at line 80 of file BinarySystem.h.
|
private |
The rate at which the orbit gains energy due to tides.
Definition at line 80 of file BinarySystem.h.
|
private |
Estimate of the error in __orbit_power due to truncating the tidal potential eccentricity expansion.
Definition at line 80 of file BinarySystem.h.
|
private |
The torque on the orbit in the coordinate system of the outermost zone of the first body.
Definition at line 119 of file BinarySystem.h.
|
private |
An estiamte of the error in ::__orbit_torque due to truncating the eccentricity series of the tidal potential.
Definition at line 119 of file BinarySystem.h.
|
private |
The current orbital angular momentum.
Definition at line 80 of file BinarySystem.h.
|
private |
The current orbital energy.
Definition at line 80 of file BinarySystem.h.
|
private |
The current semimajor axis.
Definition at line 80 of file BinarySystem.h.
|
private |
The evolution of the semimajor axis recorded by add_to_evolution() so far.
Definition at line 64 of file BinarySystem.h.
|
mutableprivate |
The current rate at which the semiamjor axis is changing.
Definition at line 111 of file BinarySystem.h.
|
private |
The rate at which the semimajor axis evolved at each.
Definition at line 64 of file BinarySystem.h.