Planetary Orbital Evolution due to Tides
Orbital evolution of two objects experiencing tides
Evolve::BinarySystem Class Reference

Describes a system of two bodies orbiting each other. More...

#include <BinarySystem.h>

+ Inheritance diagram for Evolve::BinarySystem:
+ Collaboration diagram for Evolve::BinarySystem:

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 DissipatingBodyprimary () const
 Returns the primary body in the system (const). More...
 
const DissipatingBodysecondary () 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 CombinedStoppingConditionstopping_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...
 

Detailed Description

Describes a system of two bodies orbiting each other.

The following variable are referred to throughout:

  • age: Gyr
  • a: semimajor axis in \(R_\odot\) (average distance between the bodies)
  • e: the eccentricity of the orbit
  • Eorb: the energy of the orbit (potential + kinetic) treating the two bodies as point masses in \(GM_\odot^2/R_\odot\).
  • L: angular momentum of the orbit treating the two bodies as point masses in \(M_\odot\cdot R_\odot^2 \cdot \mathrm{rad}/\mathrm{day}\)
  • \(\theta^b_i\): inclination of the i-th zone of the b-th body relative to the orbit in radians.
  • \(\omega^b_i\): argument of periapsis of the orbit for the i-th zone of the b-th body with a plane of reference perpendicular to that zone's spin.
  • \(S^b_i\): spin angular momentum of the i-th zone of the b-th body in \(M_\odot\cdot R_\odot^2 \cdot \mathrm{rad}/\mathrm{day}\).

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.

Constructor & Destructor Documentation

◆ BinarySystem()

Evolve::BinarySystem::BinarySystem ( DissipatingBody body1,
DissipatingBody body2,
const std::string &  system_name = "" 
)
inline

Construct a binary system.

Parameters
body1The first body in the system. Assumed to always be there, so for a star-planet system this should be the star.
body2The second body in the system, initially may not be there and later may be engulfed by the first body.
system_nameThe name of the system.

Definition at line 763 of file BinarySystem.h.

Member Function Documentation

◆ above_lock_fraction()

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.

Parameters
locked_zone_indexThe index within the list of locked zones of the zone we need.
entryWhat 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_indexThe index of the zone whose quantity we want the derivative w.r.t. for zone specific quantities.
secondary_radiusIf deriv==Dissipation::RADIUS this argument determines the body whose radius the derivative is w.r.t.

Definition at line 2073 of file BinarySystem.cpp.

◆ above_lock_fractions_deriv()

Eigen::VectorXd Evolve::BinarySystem::above_lock_fractions_deriv ( Dissipation::QuantityEntry  entry,
DissipatingBody body,
unsigned  zone_index 
)
private

Calculates derivatives of the above lock fractions w.r.t. quantities of non-surface zones.

Parameters
entryThe derivative to calculate. Only INCLINATION, PERIAPSIS, MOMENT_OF_INERTIA and SPIN_ANGMOM derivatives are supported.
bodyThe body whose zone's inclination we want the derivative w.r.t.
zone_indexThe zone whose inclination we are finding the derivative w.r.t.

Definition at line 700 of file BinarySystem.cpp.

◆ above_lock_problem_deriv_correction()

void Evolve::BinarySystem::above_lock_problem_deriv_correction ( Dissipation::QuantityEntry  entry,
bool  body1_deriv,
Eigen::MatrixXd &  matrix,
Eigen::VectorXd &  rhs 
) const
private

Makes corrections to the matrix and RHS to accomodate the given derivative for the linear problem that defines the above fractions.

Parameters
entryThe derivative being calculated
body1_derivFor zone-specific quantities, are we differentiating w.r.t. a zone part of body1?
matrixThe matrix to update.
rhsThe RHS vector to update.

Definition at line 501 of file BinarySystem.cpp.

◆ add_body_rate_deriv()

template<typename VALUE_TYPE >
void Evolve::BinarySystem::add_body_rate_deriv ( const DissipatingBody body,
VALUE_TYPE(DissipatingBody::*)(Dissipation::QuantityEntry, unsigned, const Eigen::VectorXd &) const  func,
std::valarray< VALUE_TYPE > &  orbit_rate_deriv,
unsigned  offset 
) const
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.

Parameters
bodyThe body whose contribution we wish to add.
funcThe quantity we are collecting. Should be either body.tidal_orbit_power or body.tidal_orbit_torque.
orbit_rate_derivThe array to add to. The indices match the parameters being evolved.
offsetIf 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.

◆ add_to_evolution()

void Evolve::BinarySystem::add_to_evolution ( )
virtual

Appends the state defined by last configure(), to the evolution.

Definition at line 2353 of file BinarySystem.cpp.

◆ age()

double Evolve::BinarySystem::age ( ) const
inline

Returns the present age of the system in Gyr.

Definition at line 833 of file BinarySystem.h.

◆ angle_evolution_age_deriv()

void Evolve::BinarySystem::angle_evolution_age_deriv ( DissipatingBody body,
unsigned  zone_ind,
double  sin_inc,
double  cos_inc,
unsigned  locked_zone_ind,
double &  inclination,
double &  periapsis 
) const
private

Computes the partial age derivatives of the inclination or periapsis evolutiono equations.

Parameters
bodyThe body whose zone to calculate the age derivative for.
zone_indThe index of the zone in body to calculate the age derivative for.
sin_incThe sin of the inclination between the zone and the orbit.
cos_incThe cos of the inclination between the zone and the orbit.
locked_zone_indThe index of the current zone in the list of locked zones (used only if zone is locked, ignored otherwise)
inclinationDestination overwritten by the age derivative of the inclination evolution rate.
periapsisDestination overwritten by the age derivative of the periapsis evolution rate.

Definition at line 1209 of file BinarySystem.cpp.

◆ angle_evolution_orbit_deriv()

void Evolve::BinarySystem::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
private

Computes the partial derivatives w.r.t. semimamjor axis or eccentricity of the inclination and periapsis evolution equations.

Parameters
entryThe derivative to compute, should be either Dissipation::SEMIMAJOR or Dissipation::ECCENTRICITY
angmom_derivThe derivative of the orbital angular momentum w.r.t. deriv.
bodyThe body whose zone's evolution we are differentiating.
zone_indThe index of the zone whose evolution we are differentiating.
sin_incThe sin of the inclination between the zone and the orbit.
cos_incThe cos of the inclination between the zone and the orbit.
locked_zone_indThe index of the current zone in the list of locked zones (used only if zone is locked, ignored otherwise)
inclinationOverwritten by the derivative of the inclination evolution rate.
periapsisOverwritten by the derivative of the periapsis evolution rate.

Definition at line 1287 of file BinarySystem.cpp.

◆ binary_differential_equations()

int Evolve::BinarySystem::binary_differential_equations ( double *  differential_equations,
bool  expansion_error 
) const
private

The differential equations for a system with both bodies present.

Also updates the current evolution rate private members if expansion_error is false.

Parameters
differential_equationsOn output is set to the rates of change of the evolution variables. See differintal_equations() for details.
expansion_errorIf 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.

◆ binary_jacobian()

void Evolve::BinarySystem::binary_jacobian ( double *  param_derivs,
double *  age_derivs 
) const
private

Calculates the jacobian for the evolution of a system of two bodies.

Parameters
param_derivsOn output is set to the Jacobian.
age_derivsOn output is set to the partial age derivatives of the evolution equations.

Definition at line 1679 of file BinarySystem.cpp.

◆ calculate_above_lock_fractions()

void Evolve::BinarySystem::calculate_above_lock_fractions ( Eigen::VectorXd &  fractions,
Dissipation::QuantityEntry  entry = Dissipation::NO_DERIV,
bool  body1_deriv = true 
)
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.

Parameters
fractionsThe vector to fill with the calculated fractions.
entryThe derivative of the above lock fractions to calculate. For zone-dependent quantities the derivatives are w.r.t. the surface zone quantity.
body1_derivOnly 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.

◆ change_e_order()

virtual void Evolve::BinarySystem::change_e_order ( unsigned  new_e_order)
inlinevirtual

Change the eccentricity expansion order for all dissipative zones.

Parameters
new_e_orderThe new eccentricity expansion order.

Definition at line 1060 of file BinarySystem.h.

◆ check_for_lock()

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.

Parameters
orbital_freq_multThe multiplier of the orbital frequency at the lock
spin_freq_multThe multiplier of the spin frequency of the potentially locked zone at the lock.
body_indexThe body whose zone is being locked.
zone_indexThe zone being unlocked.
directionThe direction in which the lock will be left if it does not hold.

Definition at line 2169 of file BinarySystem.cpp.

◆ configure() [1/2]

int Evolve::BinarySystem::configure ( bool  initialize,
double  age,
double  semimajor,
double  eccentricity,
const double *  spin_angmom,
const double *  inclination,
const double *  periapsis,
Core::EvolModeType  evolution_mode 
)
virtual

Sets the current state of the system.

Parameters
initializeIs this the first time configure() is invoked?
ageThe age to set the system to.
semimajorThe semimajor axis of the orbit.
eccentricityThe eccentricity of the orbit.
spin_angmomThe spin angular momenta of the zones of the bodies (body 1 first, outermost zone to innermost, followed by body 2).
inclinationThe inclinations of the zones of the bodies (same order as spin_angmom). The surface zone inclination must be omitted for single body systems.
periapsisThe 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_modeThe evolution mode to assume.

Definition at line 1912 of file BinarySystem.cpp.

◆ configure() [2/2]

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.

Parameters
initializeIs this the first time configure() is invoked?
ageThe current age of the system.
parametersThe parameters being evolved. See differential_equations().
evolution_modeThe evolution mode to assume for the system.

Definition at line 2005 of file BinarySystem.cpp.

◆ differential_equations()

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:

  • Calls configure().
  • Updates instantaneous evolution rates if expansion_error is false.
Parameters
ageThe current age of the system.
parametersContains 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:

  1. \(S^0_i,\ i=1\ldots\) - all zones' spin axes are assumed aligned

If evolution mode is SINGLE

  1. \(theta^0_i,\ i=1\ldots\) ( \(\theta^0_0\equiv0\))
  2. \(\omega^0_i,\ i=1\ldots\) ( \(\omega^0_0\equiv0\))
  3. \(S^0_i,\ i=0,1,\ldots\)

The TABULATION evolution mode is not allowed.

If evolution mode is BINARY

  1. a^6.5 if no zones are locked, a otherwise
  2. e
  3. \(\theta^0_i,\ i=0,\ldots\)
  4. \(\theta^1_i,\ i=0,\ldots\)
  5. \(\omega^0_i,\ i=1,\ldots\) ( \(\omega^0_0\equiv0\))
  6. \(\omega^1_i,\ i=0,\ldots\)
  7. \(S^0_i,\ i=0\ldots\) with zones in spin-orbit lock omitted
  8. \(S^1_i,\ i=0\ldots\) with zones in spin-orbit lock omitted
    Parameters
    evolution_modeThe evolution mode to assume for the system.
    differential_equationsOn 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_errorIf 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.

◆ eccentricity()

double Evolve::BinarySystem::eccentricity ( ) const
inline

The current eccentricity of the system.

Definition at line 855 of file BinarySystem.h.

◆ eccentricity_evolution() [1/2]

double Evolve::BinarySystem::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
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.

Parameters
orbit_powerSee semimajor_evolution()
orbit_angmom_gainThe 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_derivIf 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_derivIf 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_derivIf true the derivative calculated is assumed to be w.r.t. the semimajor axis.

Definition at line 430 of file BinarySystem.cpp.

◆ eccentricity_evolution() [2/2]

const std::list<double>& Evolve::BinarySystem::eccentricity_evolution ( ) const
inline

The tabulated evolution of the eccentricity so far.

Definition at line 1052 of file BinarySystem.h.

◆ eccentricity_evolution_expansion_error()

double Evolve::BinarySystem::eccentricity_evolution_expansion_error ( ) const
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.

◆ eccentricity_evolution_rate()

const std::list<double>& Evolve::BinarySystem::eccentricity_evolution_rate ( ) const
inline

The tabulated evolution of the eccentricity so far.

Definition at line 1056 of file BinarySystem.h.

◆ eccentricity_jacobian()

void Evolve::BinarySystem::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
private

Computes the row of the jacobian corresponding to the eccentricity.

Parameters
orbit_power_derivThe derivatives of the orbit energy gain w.r.t. the evolution variables and age (last entry).
orbit_angmom_gain_derivThe derivatives of the orbit angular momentum gain w.r.t. the evolution variables and age (last entry).
a6p5Is the first variable \(a^{6.5}\) instead of a?
param_derivsThe location to fill with the deriavtives of the eccentricity evolution equation w.r.t. the evolution variables.
age_derivThe location to set to the age derivative of the eccentricity evolution equation.

Definition at line 1182 of file BinarySystem.cpp.

◆ eccentricity_order()

unsigned Evolve::BinarySystem::eccentricity_order ( ) const
virtual

To what order should eccentricity expansion be performed for the given value of the eccentricity.

Definition at line 2408 of file BinarySystem.cpp.

◆ evolution_mode()

Core::EvolModeType Evolve::BinarySystem::evolution_mode ( )
inline

The evolution mode of last call to configure().

Definition at line 996 of file BinarySystem.h.

◆ fill_above_lock_fractions_deriv()

void Evolve::BinarySystem::fill_above_lock_fractions_deriv ( )
private

Fills the __above_lock_fractinos_*_deriv members.

Definition at line 779 of file BinarySystem.cpp.

◆ fill_binary_orbit()

void Evolve::BinarySystem::fill_binary_orbit ( std::valarray< double > &  orbit) const
private

Implements fill_orbit() for BINARY evolution mode.

Definition at line 1855 of file BinarySystem.cpp.

◆ fill_locked_surface_orbit()

void Evolve::BinarySystem::fill_locked_surface_orbit ( std::valarray< double > &  orbit) const
private

Implements fill_orbit() for LOCKED_SURFACE_SPIN evolution mode.

Definition at line 1839 of file BinarySystem.cpp.

◆ fill_orbit()

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.

Parameters
orbitThe orbit to fill (resized as necessary).

Definition at line 2062 of file BinarySystem.cpp.

◆ fill_orbit_angmom_gain_deriv()

void Evolve::BinarySystem::fill_orbit_angmom_gain_deriv ( std::valarray< double > &  orbit_angmom_gain_deriv) const
private

Computes the derivatives w.r.t. the evolution quantities of the orbit angular momentum gain.

Parameters
orbit_angmom_gain_derivLocation to fill.

Definition at line 1122 of file BinarySystem.cpp.

◆ fill_orbit_power_deriv()

void Evolve::BinarySystem::fill_orbit_power_deriv ( std::valarray< double > &  orbit_power_deriv) const
private

Computes the derivatives w.r.t. the evolution quantities of the orbit energy gain.

Parameters
orbit_power_derivLocation to fill.

Definition at line 1107 of file BinarySystem.cpp.

◆ fill_orbit_torque_and_power()

void Evolve::BinarySystem::fill_orbit_torque_and_power ( )
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.

◆ fill_orbit_torque_deriv()

void Evolve::BinarySystem::fill_orbit_torque_deriv ( Dissipation::QuantityEntry  entry,
DissipatingBody body,
unsigned  zone_ind,
std::valarray< Eigen::Vector3d > &  orbit_torque_deriv 
) const
private
Parameters
entryThe derivatives to compute. Sholud be one of: Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::MOMENT_OF_INERTIA, Dissipation::SPIN_ANGMOM
bodyThe body whose zone's coordinate system we are expressing the torque in.
zone_indThe index of the zone whose coordinate system we are expressing the torque in.
orbit_torque_derivOn 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.

◆ fill_single_body_jacobian()

void Evolve::BinarySystem::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
private

Fills the jacobian for a system consisting of one isolated body.

Parameters
inclination_param_derivsThe rows of the jacobian corresponding to the zone inclinations.
periapsis_param_derivsThe rows of the jacobian corresponding to the zone periapses.
angmom_param_derivsThe rows of the jacobian corresponding to the zone angular momenta.
inclination_age_derivsThe part of the age derivatives array corresponding to the zone inclinations.
periapsis_age_derivsThe part of the age derivatives array corresponding to the zone periapses.
angmom_age_derivsThe part of the age derivatives array corresponding to the zone angular momenta.

Definition at line 229 of file BinarySystem.cpp.

◆ fill_single_orbit()

void Evolve::BinarySystem::fill_single_orbit ( std::valarray< double > &  orbit) const
private

Implements fill_orbit() for SINGLE evolution mode.

Definition at line 1886 of file BinarySystem.cpp.

◆ fill_zone_torque_deriv()

void Evolve::BinarySystem::fill_zone_torque_deriv ( Dissipation::QuantityEntry  entry,
DissipatingBody body,
unsigned  zone_ind,
std::valarray< Eigen::Vector3d > &  zone_torque_deriv 
) const
private

Calculates the derivatives of the torque on a zone w.r.t. a zone-specific quantity.

Parameters
entryThe derivatives to compute. Sholud be one of: Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::MOMENT_OF_INERTIA, Dissipation::SPIN_ANGMOM
bodyThe body whose zone's evolution we are differentiating.
zone_indThe index of the zone whose evolution we are differentiating.
zone_torque_derivOverwritten 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.

◆ find_locked_zones()

void Evolve::BinarySystem::find_locked_zones ( )
private

Fills the __locked_zones list.

Definition at line 13 of file BinarySystem.cpp.

◆ get_name()

const std::string Evolve::BinarySystem::get_name ( ) const
inline

Returns the name of the system.

Definition at line 782 of file BinarySystem.h.

◆ inclination_evolution_zone_derivs()

void Evolve::BinarySystem::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
private

Calculates the derivatives of an inclination evolution equation w.r.t. a zone specific quentity.

Parameters
entryThe derivatives to compute. Sholud be one of: Dissipation::INCLINATION, Dissipation::PERIAPSIS, Dissipation::MOMENT_OF_INERTIA, Dissipation::SPIN_ANGMOM
bodyThe body whose zone's evolution we are differentiating.
zone_indThe index of the zone whose evolution we are differentiating.
zone_x_torque_aboveThe x component of the torque on the zone (above a potential spin orbit lock). Ignored if the zone is not locked.
zone_x_torque_belowThe x component of the torque on the zone (below a potential spin orbit lock).
zone_torque_derivTo be computed by fill_zone_torque_deriv().
orbit_torqueThe torque on the orbit in the coordinate system of the zone defined by the above 2 arguments.
orbit_torque_derivThe result of fill_orbit_torque_deriv.
above_frac_derivThe derivatives of the above_lock fractions.
sin_incThe sin of the inclination between the zone and the orbit.
cos_incThe cos of the inclination between the zone and the orbit.
locked_zone_indThe index of the current zone in the list of locked zones (used only if zone is locked, ignored otherwise)
resultOverwritten by the result.

Definition at line 1462 of file BinarySystem.cpp.

◆ jacobian()

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()!

Parameters
ageThe age of the system.
parametersSee differential_equations()
evolution_modeThe evolution mode to assume for the system.
param_derivsThe matrix of partial derivatives of the evolution rates for the parameters w.r.t. each of the parameters.
age_derivsThe partialderivatives of the evolution rates for the parameters w.r.t. age.

Definition at line 2146 of file BinarySystem.cpp.

◆ locked_surface_differential_equations()

int Evolve::BinarySystem::locked_surface_differential_equations ( double *  evolution_rates,
bool  expansion_error 
) const
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().

Parameters
evolution_ratesOn output is set to the rates of change of \(S^0_i\).
expansion_errorIf 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.

◆ locked_surface_jacobian()

void Evolve::BinarySystem::locked_surface_jacobian ( double *  param_derivs,
double *  age_derivs 
) const
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().

Parameters
param_derivsOn output is set to the Jacobian.
age_derivsOn output is set to the partial age derivatives of the evolution equations.

Definition at line 68 of file BinarySystem.cpp.

◆ minimum_separation()

double Evolve::BinarySystem::minimum_separation ( bool  deriv = false) const
virtual

Smallest semimajor axis at which the secondary can survive for the latest system configuration.

By default returns the larger of:

  • primary radius.
  • 2.44*(secondary radius)*((primary mass)/(secondary mass))^(1/3)
Parameters
derivIf true the rate of change (per Gyr) is returned.

Definition at line 2248 of file BinarySystem.cpp.

◆ next_stop_age()

double Evolve::BinarySystem::next_stop_age ( ) const
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.

◆ number_locked_zones()

unsigned Evolve::BinarySystem::number_locked_zones ( ) const
inline

How many zones on either body are currently locked.

Definition at line 848 of file BinarySystem.h.

◆ number_zones()

unsigned Evolve::BinarySystem::number_zones ( ) const
inline

The total number of zones in both system bodies.

Definition at line 842 of file BinarySystem.h.

◆ periapsis_evolution_zone_derivs()

void Evolve::BinarySystem::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
private

Computes the derivatives of the periapsis evolution w.r.t. a zone specific quantity for all zones.

Parameters
entryThe derivative to compute.
bodyThe body whose zone's periapsis evolution to differentiate.
zone_indThe zone whose periapsis evolution to differentiate.
zone_y_torque_aboveThe y component of the torque on the zone above a lock if locked.
zone_y_torque_belowThe y component of the torque on the zone (assuming below a lock if locked).
zone_torque_derivThe derivative of the zone torque.
orbit_y_torqueThe y torque on the orbit.
orbit_torque_derivThe derivative of the orbit torque.
above_frac_derivThe derivatives of the above lock fractions
sin_incThe sin of the inclination between the zone and the orbit.
cos_incThe cos of the inclination between the zone and the orbit.
locked_zone_indThe index in the list of locked zones of this zone if it is locked.
resultOverwritten by the result. Should have the correct size.

Definition at line 1542 of file BinarySystem.cpp.

◆ primary()

const DissipatingBody& Evolve::BinarySystem::primary ( ) const
inline

Returns the primary body in the system (const).

Definition at line 836 of file BinarySystem.h.

◆ reached_critical_age()

void Evolve::BinarySystem::reached_critical_age ( double  age)
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.

◆ release_lock()

void Evolve::BinarySystem::release_lock ( unsigned  locked_zone_index,
short  direction 
)
virtual

Releases the lock to one of the locked zones.

Parameters
locked_zone_indexThe index within the list of locked zones of the zone to unlock.
directionThe direction in which the zone's spin will evolve in the future relative to the lock.

Definition at line 2321 of file BinarySystem.cpp.

◆ reset_evolution()

void Evolve::BinarySystem::reset_evolution ( )
virtual

Resets the evolution of the system.

Definition at line 2363 of file BinarySystem.cpp.

◆ rewind_evolution()

void Evolve::BinarySystem::rewind_evolution ( unsigned  nsteps)
virtual

Discards the last steps from the evolution.

Parameters
nstepsHow many steps of evolution to discard.

Definition at line 2373 of file BinarySystem.cpp.

◆ secondary()

const DissipatingBody& Evolve::BinarySystem::secondary ( ) const
inline

Returns the secondary body in the system (const).

Definition at line 839 of file BinarySystem.h.

◆ secondary_died()

void Evolve::BinarySystem::secondary_died ( )
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.

◆ semimajor()

double Evolve::BinarySystem::semimajor ( ) const
inline

The current semimajor axis of the system.

Definition at line 852 of file BinarySystem.h.

◆ semimajor_evolution() [1/2]

double Evolve::BinarySystem::semimajor_evolution ( double  orbit_power,
double  orbit_power_deriv = Core::NaN 
) const
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.

Parameters
orbit_powerThe 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_derivIf 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.

◆ semimajor_evolution() [2/2]

const std::list<double>& Evolve::BinarySystem::semimajor_evolution ( ) const
inline

The tabulated evolution of the semimajor axis so far.

Definition at line 1044 of file BinarySystem.h.

◆ semimajor_evolution_expansion_error()

double Evolve::BinarySystem::semimajor_evolution_expansion_error ( ) const
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.

◆ semimajor_evolution_rate()

const std::list<double>& Evolve::BinarySystem::semimajor_evolution_rate ( ) const
inline

The tabulated evolution of the semimajor axis so far.

Definition at line 1048 of file BinarySystem.h.

◆ semimajor_jacobian()

void Evolve::BinarySystem::semimajor_jacobian ( const std::valarray< double > &  orbit_power_deriv,
bool  a6p5,
double *  param_derivs,
double &  age_deriv 
) const
private

Computes the row of the jacobian corresponding to the semimajor axis.

Parameters
orbit_power_derivThe derivatives of the orbit energy gain w.r.t. the evolution variables and age (last entry).
a6p5Is the first variable \(a^{6.5}\) instead of a?
param_derivsThe location to fill with the deriavtives of the semimajor evolution equation w.r.t. the evolution variables.
age_derivThe location to set to the age derivative of the semimajor evolution equation.

Definition at line 1165 of file BinarySystem.cpp.

◆ single_body_differential_equations()

int Evolve::BinarySystem::single_body_differential_equations ( double *  evolution_rates,
bool  expansion_error 
) const
private

Differential equations for the rotation of the zones of body 0 if no other body is present.

configure() must already have been called.

Parameters
evolution_ratesOn outputs is set to the rate of change of the orbital parameters.
expansion_errorIf 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.

◆ single_body_jacobian()

void Evolve::BinarySystem::single_body_jacobian ( double *  param_derivs,
double *  age_derivs 
) const
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().

Parameters
param_derivsOn output is set to the Jacobian.
age_derivsOn output is set to the partial age derivatives of the evolution equations.

Definition at line 394 of file BinarySystem.cpp.

◆ spin_angmom_evolution_zone_derivs()

void Evolve::BinarySystem::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
private
Parameters
entryThe derivative to compute.
bodyThe body whose zone's periapsis evolution to differentiate.
zone_indThe zone whose periapsis evolution to differentiate.
zone_z_torque_aboveThe z component of the torque on the zone above a lock if locked.
zone_z_torque_belowThe z component of the torque on the zone (assuming below a lock if locked).
zone_torque_derivThe derivative of the zone torque.
above_frac_derivThe derivatives of the above lock fractions
locked_zone_indThe index in the list of locked zones of this zone if it is locked.
resultOverwritten by the result. Should have the correct size.

Definition at line 1635 of file BinarySystem.cpp.

◆ stopping_conditions()

CombinedStoppingCondition * Evolve::BinarySystem::stopping_conditions ( )
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.

◆ update_above_lock_fractions()

void Evolve::BinarySystem::update_above_lock_fractions ( )
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.

Member Data Documentation

◆ __above_lock_fractions

std::valarray<Eigen::VectorXd> Evolve::BinarySystem::__above_lock_fractions
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.

◆ __above_lock_fractions_angmom_deriv

std::valarray<Eigen::VectorXd> Evolve::BinarySystem::__above_lock_fractions_angmom_deriv
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.

◆ __above_lock_fractions_body2_radius_deriv

Eigen::VectorXd Evolve::BinarySystem::__above_lock_fractions_body2_radius_deriv
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.

◆ __above_lock_fractions_decomp

Eigen::ColPivHouseholderQR<Eigen::MatrixXd> Evolve::BinarySystem::__above_lock_fractions_decomp
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.

◆ __above_lock_fractions_inclination_deriv

std::valarray<Eigen::VectorXd> Evolve::BinarySystem::__above_lock_fractions_inclination_deriv
private

The derivatives of the above lock fractions w.r.t. the inclinations of the zones.

Definition at line 136 of file BinarySystem.h.

◆ __above_lock_fractions_inertia_deriv

std::valarray<Eigen::VectorXd> Evolve::BinarySystem::__above_lock_fractions_inertia_deriv
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.

◆ __above_lock_fractions_periapsis_deriv

std::valarray<Eigen::VectorXd> Evolve::BinarySystem::__above_lock_fractions_periapsis_deriv
private

The derivatives of the above lock fractions w.r.t. the periapses of the zones.

Definition at line 136 of file BinarySystem.h.

◆ __age

double Evolve::BinarySystem::__age
private

The present age of the stellar system in Gyrs.

Definition at line 80 of file BinarySystem.h.

◆ __body1

DissipatingBody& Evolve::BinarySystem::__body1
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.

◆ __body2

DissipatingBody & Evolve::BinarySystem::__body2
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.

◆ __eccentricity

double Evolve::BinarySystem::__eccentricity
private

The current eccentricity.

Definition at line 80 of file BinarySystem.h.

◆ __eccentricity_evolution

std::list<double> Evolve::BinarySystem::__eccentricity_evolution
private

The evolution of the eccentricity recorded by add_to_evolution() so far.

Definition at line 64 of file BinarySystem.h.

◆ __eccentricity_rate

double Evolve::BinarySystem::__eccentricity_rate
mutableprivate

The current rate at which the eccentricity is changing.

Definition at line 111 of file BinarySystem.h.

◆ __eccentricity_rate_evolution

std::list<double> Evolve::BinarySystem::__eccentricity_rate_evolution
private

The rate at which the eccentricity evolved at each.

Definition at line 64 of file BinarySystem.h.

◆ __evolution_mode

Core::EvolModeType Evolve::BinarySystem::__evolution_mode
private

The evolution mode from the last call to configure();.

Definition at line 126 of file BinarySystem.h.

◆ __locked_zones

std::list<unsigned> Evolve::BinarySystem::__locked_zones
private

A list of indices of locked zones.

Definition at line 129 of file BinarySystem.h.

◆ __name

std::string Evolve::BinarySystem::__name
private

The name of the binary system (e.g. "HAT-P-20")

Definition at line 59 of file BinarySystem.h.

◆ __orbit_angmom_gain

double Evolve::BinarySystem::__orbit_angmom_gain
private

The rate at which the orbit gains angular momentum due to tides.

Definition at line 80 of file BinarySystem.h.

◆ __orbit_angmom_gain_expansion_error

double Evolve::BinarySystem::__orbit_angmom_gain_expansion_error
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.

◆ __orbit_power

double Evolve::BinarySystem::__orbit_power
private

The rate at which the orbit gains energy due to tides.

Definition at line 80 of file BinarySystem.h.

◆ __orbit_power_expansion_error

double Evolve::BinarySystem::__orbit_power_expansion_error
private

Estimate of the error in __orbit_power due to truncating the tidal potential eccentricity expansion.

Definition at line 80 of file BinarySystem.h.

◆ __orbit_torque

Eigen::Vector3d Evolve::BinarySystem::__orbit_torque
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.

◆ __orbit_torque_expansion_error

Eigen::Vector3d Evolve::BinarySystem::__orbit_torque_expansion_error
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.

◆ __orbital_angmom

double Evolve::BinarySystem::__orbital_angmom
private

The current orbital angular momentum.

Definition at line 80 of file BinarySystem.h.

◆ __orbital_energy

double Evolve::BinarySystem::__orbital_energy
private

The current orbital energy.

Definition at line 80 of file BinarySystem.h.

◆ __semimajor

double Evolve::BinarySystem::__semimajor
private

The current semimajor axis.

Definition at line 80 of file BinarySystem.h.

◆ __semimajor_evolution

std::list<double> Evolve::BinarySystem::__semimajor_evolution
private

The evolution of the semimajor axis recorded by add_to_evolution() so far.

Definition at line 64 of file BinarySystem.h.

◆ __semimajor_rate

double Evolve::BinarySystem::__semimajor_rate
mutableprivate

The current rate at which the semiamjor axis is changing.

Definition at line 111 of file BinarySystem.h.

◆ __semimajor_rate_evolution

std::list<double> Evolve::BinarySystem::__semimajor_rate_evolution
private

The rate at which the semimajor axis evolved at each.

Definition at line 64 of file BinarySystem.h.


The documentation for this class was generated from the following files: