A base class for any body contributing to tidal dissipation. More...
#include <DissipatingBody.h>
Public Member Functions | |
DissipatingBody () | |
Some initializations for new objects. More... | |
virtual void | configure (bool initialize, double age, double companion_mass, double semimajor, double eccentricity, const double *spin_angmom, const double *inclination=NULL, const double *periapsis=NULL, bool locked_surface=false, bool zero_outer_inclination=false, bool zero_outer_periapsis=false) |
Defines the orbit this body is in. More... | |
void | lock_zone_spin (unsigned zone_index, int orbital_frequency_multiplier, int spin_frequency_multiplier) |
void | unlock_zone_spin (unsigned zone_index, short direction) |
Releases the given zone from a spin-orbit lock. More... | |
unsigned | number_locked_zones () const |
The number of zones currently in a spin-orbit lock. More... | |
Eigen::Vector3d | nontidal_torque (unsigned zone_index, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, int deriv_zone=0) const |
External torque acting on a single zone (last calculate_torques_power()). More... | |
const Eigen::Vector3d & | tidal_torque (unsigned zone_index, bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const |
Tidal torque acting on the given zone (last calculate_torques_power()). More... | |
double | tidal_power (unsigned zone_index, bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const |
Tidal power dissipated in the given zone. More... | |
void | set_above_lock_fractions (std::valarray< Eigen::VectorXd > &above_lock_fractions) |
Corrects the tidal orbit energy gain and angular momentum gain for locked zones. More... | |
double | tidal_orbit_power (Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, unsigned deriv_zone_index=0, const Eigen::VectorXd &above_lock_fraction_deriv=Eigen::VectorXd()) const |
Rate of increase of the orbital energy due to tides in this body (last calculate_torques_power()). More... | |
Eigen::Vector3d | tidal_orbit_torque (Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, unsigned deriv_zone_index=0, const Eigen::VectorXd &above_lock_fraction_deriv=Eigen::VectorXd()) const |
The torque on the orbit due to tidal dissipation in the body. More... | |
Eigen::Vector3d | tidal_orbit_torque (const DissipatingZone &reference_zone, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, unsigned deriv_zone_index=0, const Eigen::VectorXd &above_lock_fraction_deriv=Eigen::VectorXd()) const |
Same as tidal_orbit_torque(Dissipation::QuantityEntry, unsigned, const Eigen::VectorXd &) but allow specifying the zone whose coordinate system to use. More... | |
virtual unsigned | number_zones () const =0 |
The number of zones the body consists of. More... | |
virtual const DissipatingZone & | zone (unsigned zone_index) const =0 |
A modifiable reference to one of the body's zones. More... | |
virtual DissipatingZone & | zone (unsigned zone_index)=0 |
A modifiable reference to one of the body's zones. More... | |
virtual Eigen::Vector3d | angular_momentum_coupling (unsigned top_zone_index, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, bool with_respect_to_top=false) const =0 |
Coupling torque for two neighboring zones in the coordinate system of the top zone. More... | |
virtual double | angular_momentum_loss (Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV) const =0 |
Rate of angular momentum loss by the top zone of the body and its derivatives. More... | |
double | radius (int deriv_order=0) const |
The current radius or its derivative with age of the body. More... | |
double | mass () const |
The mass of the body (constant with age). More... | |
double | spin_frequency () const |
The surface spin freuqency of the body. More... | |
double | surface_lock_frequency () const |
Angular velocity of the surface zone when locked (assumed constant). More... | |
void | set_surface_lock_frequency (double frequency) |
Sets the frequency at which the surface is locked (if any). More... | |
virtual void | add_to_evolution () |
Appends the state defined by last configure(), to the evolution. More... | |
virtual void | rewind_evolution (unsigned nsteps) |
Discards the last steps from the evolution. More... | |
virtual void | reset_evolution () |
Discards all evolution. More... | |
virtual CombinedStoppingCondition * | stopping_conditions (BinarySystem &system, bool primary) |
Conditions detecting the next possible discontinuities in the evolution due to this body. More... | |
virtual void | spin_jumped () |
Notifies the body that its spin just discontinously jumped. More... | |
virtual void | reached_critical_age (double) |
Change the body 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 change in one of the bodies. More... | |
virtual void | change_e_order (unsigned new_e_order, BinarySystem &system, bool primary) |
Change the eccentricity expansion order for all dissipative zones. More... | |
virtual | ~DissipatingBody () |
Virtual destructor. More... | |
Private Member Functions | |
double | normalize_torques (double companion_mass, double semimajor, double orbital_frequency) |
Scales the dimensionless torques as appropriate and corrects the relevant derivatives returning the normalization used. More... | |
void | collect_orbit_rates (double orbital_frequency, double torque_norm) |
Calcuates the total energy and angular momentum loss rates for the orbit due to this body's tidal dissipation. More... | |
void | calculate_orbit_rate_corrections () |
Calculates the correction the orbital energy gain and torque due to switching locked zones from fully below to fully above. More... | |
void | angular_momentum_transfer (const DissipatingZone &outer_zone, const DissipatingZone &inner_zone, Eigen::Vector3d &outer_angmom_gain, Eigen::Vector3d &inner_angmom_gain, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, bool with_respect_to_outer=false) const |
Angular momentum transfer between two neighboring zones due to moving zone boundaries. More... | |
Eigen::Vector3d | angular_momentum_transfer_from_top (unsigned zone_index, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, bool with_respect_to_outer=false) const |
Rate of angular momentum transfer (or its derivatives) to a zone due to its top boundary moving. More... | |
Eigen::Vector3d | angular_momentum_transfer_from_bottom (unsigned zone_index, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, bool with_respect_to_inner=false) const |
Rate of angular momentum transfer (or its derivatives) to a zone due to its bottom boundary moving. More... | |
Eigen::Vector3d | angular_momentum_transfer_to_zone (unsigned zone_index, Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV, int deriv_zone=0) const |
Rate of angular momentum transfer (or its derivatives) to a zone due to moving zone boundaries. More... | |
void | calculate_nontidal_torques () |
Calculates the non-tidal torques on all zones. More... | |
void | correct_orbit_power (Eigen::VectorXd &above_lock_fractions_age_deriv, Eigen::VectorXd &above_lock_fractions_semimajor_deriv, Eigen::VectorXd &above_lock_fractions_eccentricity_deriv, Eigen::VectorXd &above_lock_fractions_radius_deriv) |
Corrects __orbit_power for zone locks. More... | |
void | correct_orbit_torque (std::valarray< Eigen::VectorXd > &above_lock_fractions) |
Corrects __orbit_torque for zone locks. More... | |
Private Attributes | |
double | __power_norm |
The coefficient used to normalize tidal power. More... | |
double | __orbital_frequency |
The mean angular velocity with which the orbit is traversed. More... | |
double | __dorbital_frequency_da |
The derivative of the orbital frequency w.r.t. semimajor axis. More... | |
double | __surface_lock_frequency |
The frequency at which the surface is locked (if any). More... | |
std::valarray< std::valarray< Eigen::Vector3d > > | __angular_momentum_transfer |
The rate of angular momentum transfer between two neighboring zones due to zone boundaries moving. More... | |
std::valarray< std::valarray< Eigen::Vector3d > > | __tidal_torques_above |
Tidal torques and their derivatives on each zone for spin frequency approaching potential lock from above. More... | |
std::valarray< std::valarray< Eigen::Vector3d > > | __tidal_torques_below |
Tidal torques and their derivatives on each zone for spin frequency approaching potential lock from below. More... | |
std::vector< Dissipation::QuantityEntry > | __orbit_entries |
The quantities w.r.t. which derivatives of the orbit energy gain and torque are pre-calculated or error, if the value is __expansion_error_id. More... | |
std::valarray< double > | __orbit_power |
Total tidal power (and derivatives and error) gained by the orbit from this body. More... | |
std::valarray< double > | __orbit_power_correction |
Corrections to __orbit_power_below (undifferentiated) if single zones switch to above. More... | |
std::vector< Eigen::Vector3d > | __orbit_torque |
Torque on the orbit (and its derivatives and error) due to tides on this body. More... | |
std::vector< Eigen::Vector3d > | __orbit_torque_correction |
Corrections to __orbit_torque_below (undifferentiated) if single zones switch to above. More... | |
unsigned | __num_locked_zones |
The number of zones currently in a spin-orbit lock. More... | |
std::valarray< Eigen::VectorXd > | __above_lock_fractions |
The fractional contribution of the above the lock rates for locked zones and their derivatives. More... | |
A base class for any body contributing to tidal dissipation.
All torques are in units of \(M_\odot R_\odot^2 \mathrm{rad}\mathrm{day}^{-1}\mathrm{Gyr}^{-1}\) and tidal power is in units of \(M_\odot R_\odot^2 \mathrm{rad}^2\mathrm{day}^{-2}\mathrm{Gyr}^{-1}\)
Age derivatives are per Gyr
Masses and radii are in solar units.
Angular velocities are in rad/day.
Definition at line 39 of file DissipatingBody.h.
Evolve::DissipatingBody::DissipatingBody | ( | ) |
Some initializations for new objects.
Definition at line 422 of file DissipatingBody.cpp.
|
inlinevirtual |
Virtual destructor.
Definition at line 601 of file DissipatingBody.h.
|
virtual |
Appends the state defined by last configure(), to the evolution.
Reimplemented in Star::SaturatingSkumanichWindBody.
Definition at line 1028 of file DissipatingBody.cpp.
|
pure virtual |
Coupling torque for two neighboring zones in the coordinate system of the top zone.
Units: \(M_\odot R_\odot^2\mathrm{rad}\mathrm{day}^{-1}mathrm{Gyr}^{-1}\).
top_zone_index | The index of the outer of the two zones whose coupling is needed. |
deriv | The derivative of the coupling torque required. |
with_respect_to_top | For derivatives with respect to zone specific quantities, this determines which zone's quantity to differentiate with respect to (top zone if true, bottom zone if false). |
Implemented in Evolve::SingleTidalTermBody.
|
pure virtual |
Rate of angular momentum loss by the top zone of the body and its derivatives.
The spin frequency of the topmost zone of the body must already be set.
deriv | The derivative of the angular momentum loss required. |
|
private |
Angular momentum transfer between two neighboring zones due to moving zone boundaries.
The argument of periapsis and the inclination must be defined for both zones involved.
outer_zone | The outer zone. |
inner_zone | The inner zone. |
outer_angmom_gain | The angular momentum change for the outer zone in the outer zone's coordinate system. |
inner_angmom_gain | The angular momentum change for the inner zone in the inner zone's coordinate system. |
deriv | Derivatives with respect to inclination and periapsis can be computed, in addition to the actual transfer. It is an error to request another derivative. |
with_respect_to_outer | If deriv is not NO_DERIV, derivatives can be computed with respect to quantities of the outer zone (if this argument is true) or the inner zone (if false). |
Definition at line 187 of file DissipatingBody.cpp.
|
private |
Rate of angular momentum transfer (or its derivatives) to a zone due to its bottom boundary moving.
The result is in the coordinate system of the zone.
The configure() method must already have been called.
zone_index | The index of the zone to calculate angular momentum gain for. Must not be less than number_zones()-1. |
deriv | Whether to return the quantity or one of its derivatives. |
with_respect_to_inner | If deriv is a zone specific quantity this argument determines if derivative with respect to the quantity of the zone below (true) or this zone (false) should be calculated. |
Definition at line 306 of file DissipatingBody.cpp.
|
private |
Rate of angular momentum transfer (or its derivatives) to a zone due to its top boundary moving.
The result is in the coordinate system of the zone.
The set_orbit() method must already have been called and all zones must have their inclinations and arguments of periapsis set.
zone_index | The index of the zone to calculate angular momentum gain for. Must not be zero. |
deriv | Whether to return the quantity or one of its derivatives. |
with_respect_to_outer | If deriv is a zone specific quantity this argument determines if derivative with respect to the quantity of the zone above (true) or this zone (false) should be calculated. |
Definition at line 240 of file DissipatingBody.cpp.
|
private |
Rate of angular momentum transfer (or its derivatives) to a zone due to moving zone boundaries.
The result is in the coordinate system of the zone.
The configure() method must already have been called.
zone_index | The index of the zone to calculate angular momentum gain for. |
deriv | Whether to return the quantity or one of its derivatives. |
deriv_zone | See matching argument of external_torque() for description. |
Definition at line 372 of file DissipatingBody.cpp.
|
private |
Calculates the non-tidal torques on all zones.
|
private |
Calculates the correction the orbital energy gain and torque due to switching locked zones from fully below to fully above.
Definition at line 164 of file DissipatingBody.cpp.
|
virtual |
Change the eccentricity expansion order for all dissipative zones.
new_e_order | The new eccentricity expansion order. |
system | The system being evolved. |
primary | Is the body the primary in the system. |
Definition at line 1059 of file DissipatingBody.cpp.
|
private |
Calcuates the total energy and angular momentum loss rates for the orbit due to this body's tidal dissipation.
orbital_frequency | The orbital frequency. |
torque_norm | The normalization constant for the torques, as returned by normalize_torques(). |
Definition at line 79 of file DissipatingBody.cpp.
|
virtual |
Defines the orbit this body is in.
The inclinations and arguments of periapsis must be already set for all zones.
initialize | Is this the first time the body is configure() -ed? |
age | The age to set the body to. |
companion_mass | The mass of the second body in the system. |
semimajor | The semimajor axis of the orbit in \(R_\odot\). |
eccentricity | The eccentricity of the orbit |
spin_angmom | The spin angular momenta of the non-locked zones of the body (outermost zone to innermost). |
inclination | The inclinations of the zones of the body (same order as spin_angmom). If NULL, all inclinations are assumed zero. |
periapsis | The arguments of periapsis of the zones of the bodies (same order as spin_angmom). If NULL, all periapses are assumed zero. |
locked_surface | If true, the outermost zone's spin is assumed locked to a disk and spin_angmom is assumed to start from the next zone. |
zero_outer_inclination | If true, the outermost zone's inclination is assumed to be zero and the inclination argument is assumed to start from the next zone. |
zero_outer_periapsis | If true, the outermost zone's periapsis is assumed to be zero and the inclination argument is assumed to start from the next zone. |
Reimplemented in Star::ExponentialDecayDiffRotBody.
Definition at line 603 of file DissipatingBody.cpp.
|
private |
Corrects __orbit_power for zone locks.
above_lock_fractions_age_deriv | The derivative w.r.t. age of the above lock fractions. |
above_lock_fractions_semimajor_deriv | The derivative w.r.t. semimajor axis of the above lock fractions. |
above_lock_fractions_eccentricity_deriv | The derivative w.r.t. eccentricity of the above lock fractions. |
above_lock_fractions_radius_deriv | The derivative w.r.t. the body radius of the above lock fractions. |
Definition at line 438 of file DissipatingBody.cpp.
|
private |
Corrects __orbit_torque for zone locks.
above_lock_fractions | The same as the argument of set_above_lock_fractions. |
Definition at line 523 of file DissipatingBody.cpp.
|
inline |
Locks the given zone's spin to the orbit with the given frequency ratio.
Definition at line 314 of file DissipatingBody.h.
|
inline |
The mass of the body (constant with age).
Definition at line 534 of file DissipatingBody.h.
|
inlinevirtual |
The next age when the evolution needs to be stopped for a change in one of the bodies.
Reimplemented in Star::InterpolatedEvolutionStar.
Definition at line 586 of file DissipatingBody.h.
Eigen::Vector3d Evolve::DissipatingBody::nontidal_torque | ( | unsigned | zone_index, |
Dissipation::QuantityEntry | deriv = Dissipation::NO_DERIV , |
||
int | deriv_zone = 0 |
||
) | const |
External torque acting on a single zone (last calculate_torques_power()).
This includes torques coupling it to other zones (due to differential rotation and/or due to mass transfer) and wind loss for the surface zone.
For each zone the torques are in a coordinate system with \(\hat{x}\) along the ascending node of the orbit in the zone's equatorial plane, and \(\hat{z}\) along the zone's angular momentum.
zone_index | The index of the zone whose torque is needed. |
deriv | Which entry to return for the torque. |
deriv_zone | Since external torques depend on neighboring zones, this parameter is used to distinguish those (it is ignored if : deriv is Dissipation::NO_DERIV or a quantity which is not zone specific):
|
Definition at line 742 of file DissipatingBody.cpp.
|
private |
Scales the dimensionless torques as appropriate and corrects the relevant derivatives returning the normalization used.
companion_mass | The mass of the other object in the system. |
semimajor | The orbital semimajor axis. |
orbital_frequency | The orbital frequency (passed to avoid duplicate calculation). |
Definition at line 6 of file DissipatingBody.cpp.
|
inline |
The number of zones currently in a spin-orbit lock.
Definition at line 339 of file DissipatingBody.h.
|
pure virtual |
The number of zones the body consists of.
Implemented in Star::InterpolatedEvolutionStar, Evolve::SingleTidalTermBody, and Planet::Planet.
|
inline |
The current radius or its derivative with age of the body.
deriv_order | The order of the derivative to return. |
Definition at line 527 of file DissipatingBody.h.
|
inlinevirtual |
Change the body as necessary at the given age.
Handles things like interpolation discontinuities.
Reimplemented in Star::InterpolatedEvolutionStar, and Planet::Planet.
Definition at line 582 of file DissipatingBody.h.
|
virtual |
Discards all evolution.
Reimplemented in Star::SaturatingSkumanichWindBody.
Definition at line 1040 of file DissipatingBody.cpp.
|
virtual |
Discards the last steps from the evolution.
nsteps | How many steps of evolution to discard. |
Reimplemented in Star::SaturatingSkumanichWindBody.
Definition at line 1034 of file DissipatingBody.cpp.
void Evolve::DissipatingBody::set_above_lock_fractions | ( | std::valarray< Eigen::VectorXd > & | above_lock_fractions | ) |
Corrects the tidal orbit energy gain and angular momentum gain for locked zones.
above_lock_fractions | The fractional contributions of the above the lock rates for each zone (the vector) and their derivatives (the outer index). |
Definition at line 827 of file DissipatingBody.cpp.
|
inline |
Sets the frequency at which the surface is locked (if any).
Definition at line 549 of file DissipatingBody.h.
|
inline |
The surface spin freuqency of the body.
Definition at line 538 of file DissipatingBody.h.
|
inlinevirtual |
Notifies the body that its spin just discontinously jumped.
Reimplemented in Star::SaturatingSkumanichWindBody.
Definition at line 577 of file DissipatingBody.h.
|
virtual |
Conditions detecting the next possible discontinuities in the evolution due to this body.
Must be deleted when no longer necessary.
system | The system being evolved. |
primary | Is the body the primary in the system. |
Reimplemented in Star::SaturatingSkumanichWindBody.
Definition at line 1046 of file DissipatingBody.cpp.
|
inline |
Angular velocity of the surface zone when locked (assumed constant).
For example this could be the frequency of the stellar convective zone when locked to the disk.
Definition at line 545 of file DissipatingBody.h.
double Evolve::DissipatingBody::tidal_orbit_power | ( | Dissipation::QuantityEntry | entry = Dissipation::NO_DERIV , |
unsigned | deriv_zone_index = 0 , |
||
const Eigen::VectorXd & | above_lock_fraction_deriv = Eigen::VectorXd() |
||
) | const |
Rate of increase of the orbital energy due to tides in this body (last calculate_torques_power()).
If set_above_lock_fractions() has not been called since cofigure(), returns the rate assuming all locked zones are fully below the lock. If it has, returns the actual rate.
entry | W.r.t. that quantity is the required derivative. |
deriv_zone_index | If deriv is a zone-specific quantity, this argument determines which zone to differentiate w.r.t. |
above_lock_fraction_deriv | The derivatives of all lock fractions w.r.t. to the quantity identified by the above arguments. |
Definition at line 849 of file DissipatingBody.cpp.
Eigen::Vector3d Evolve::DissipatingBody::tidal_orbit_torque | ( | Dissipation::QuantityEntry | deriv = Dissipation::NO_DERIV , |
unsigned | deriv_zone_index = 0 , |
||
const Eigen::VectorXd & | above_lock_fraction_deriv = Eigen::VectorXd() |
||
) | const |
The torque on the orbit due to tidal dissipation in the body.
If set_above_lock_fractions() has not been called since cofigure(), returns the torue assuming all locked zones are fully below the lock. If it has, returns the actual rate.
deriv | The derivative of the angular momentum rate required. |
deriv_zone_index | If deriv is a zone-specific quantity, this argument determines which zone to differentiate w.r.t. |
above_lock_fraction_deriv | The derivatives of all lock fractions w.r.t. to the quantity identified by the above arguments. |
Definition at line 906 of file DissipatingBody.cpp.
Eigen::Vector3d Evolve::DissipatingBody::tidal_orbit_torque | ( | const DissipatingZone & | reference_zone, |
Dissipation::QuantityEntry | deriv = Dissipation::NO_DERIV , |
||
unsigned | deriv_zone_index = 0 , |
||
const Eigen::VectorXd & | above_lock_fraction_deriv = Eigen::VectorXd() |
||
) | const |
Same as tidal_orbit_torque(Dissipation::QuantityEntry, unsigned, const Eigen::VectorXd &) but allow specifying the zone whose coordinate system to use.
reference_zone | The zone whose coordinate system to express the result. |
deriv | The derivative of the angular momentum rate required. |
deriv_zone_index | If deriv is a zone-specific quantity, this argument determines which zone to differentiate w.r.t. |
above_lock_fraction_deriv | The derivatives of all lock fractions w.r.t. to the quantity identified by the above arguments. |
Definition at line 985 of file DissipatingBody.cpp.
double Evolve::DissipatingBody::tidal_power | ( | unsigned | zone_index, |
bool | above, | ||
Dissipation::QuantityEntry | entry = Dissipation::NO_DERIV |
||
) | const |
Tidal power dissipated in the given zone.
zone_index | The index of the zone whose tidal torque is required. |
above | If a zone is currently in a lock, this argument decides whetherthe torque for spin frequency above (true) or below (false) of the lock should be returned. |
entry | Which derivative of the tidal power is required. |
Definition at line 797 of file DissipatingBody.cpp.
|
inline |
Tidal torque acting on the given zone (last calculate_torques_power()).
The coordinate system is the same as for nontidal_torque().
zone_index | The index of the zone whose tidal torque is required. |
above | If a zone is currently in a lock, this argument decides whetherthe torque for spin frequency above (true) or below (false) of the lock should be returned. |
entry | Which derivative of the tidal torque is required. |
Definition at line 376 of file DissipatingBody.h.
|
inline |
Releases the given zone from a spin-orbit lock.
zone_index | The zone whose spin orbit lock to release. |
direction | The direction in which the spin will evolve in the future relative to the lock. |
Definition at line 324 of file DissipatingBody.h.
|
pure virtual |
A modifiable reference to one of the body's zones.
zone_index | The index of the zone within the body. Sequential zones are ordered from outside to inside (0 is the surface zone, number_zones()-1 is the core). |
Implemented in Star::InterpolatedEvolutionStar, Evolve::SingleTidalTermBody, and Planet::Planet.
|
pure virtual |
A modifiable reference to one of the body's zones.
zone_index | The index of the zone within the body. Sequential zones are ordered from outside to inside (0 is the surface zone, number_zones()-1 is the core). |
Implemented in Star::InterpolatedEvolutionStar, Evolve::SingleTidalTermBody, and Planet::Planet.
|
private |
The fractional contribution of the above the lock rates for locked zones and their derivatives.
Definition at line 115 of file DissipatingBody.h.
|
private |
The rate of angular momentum transfer between two neighboring zones due to zone boundaries moving.
The outer index is the index of the outer zone and each entry consists of two sub-entries giving the rate at which each of the two zones (outer first, inner second) involved gains angular momentum due to the moving boundary in each zone's coordinate system.
Definition at line 63 of file DissipatingBody.h.
|
private |
The derivative of the orbital frequency w.r.t. semimajor axis.
Definition at line 43 of file DissipatingBody.h.
|
private |
The number of zones currently in a spin-orbit lock.
Definition at line 111 of file DissipatingBody.h.
|
private |
The quantities w.r.t. which derivatives of the orbit energy gain and torque are pre-calculated or error, if the value is __expansion_error_id.
Definition at line 80 of file DissipatingBody.h.
|
private |
Total tidal power (and derivatives and error) gained by the orbit from this body.
The derivatives are only w.r.t. the non-zone specific quantities listed in __orbit_entries.
Definition at line 88 of file DissipatingBody.h.
|
private |
Corrections to __orbit_power_below (undifferentiated) if single zones switch to above.
Definition at line 88 of file DissipatingBody.h.
|
private |
Torque on the orbit (and its derivatives and error) due to tides on this body.
In the coordinate system of the topmost zone.
The derivatives are only w.r.t. the non-zone specific quantities listed in __orbit_entries.
Definition at line 101 of file DissipatingBody.h.
|
private |
Corrections to __orbit_torque_below (undifferentiated) if single zones switch to above.
It has already been transformed to the coordinate system of the surface zone.
Definition at line 101 of file DissipatingBody.h.
|
private |
The mean angular velocity with which the orbit is traversed.
Definition at line 43 of file DissipatingBody.h.
|
private |
The coefficient used to normalize tidal power.
Definition at line 43 of file DissipatingBody.h.
|
private |
The frequency at which the surface is locked (if any).
Definition at line 43 of file DissipatingBody.h.
|
private |
Tidal torques and their derivatives on each zone for spin frequency approaching potential lock from above.
See tidal_torque() for more details.
Definition at line 63 of file DissipatingBody.h.
|
private |
Tidal torques and their derivatives on each zone for spin frequency approaching potential lock from below.
See tidal_torque() for more details.
Definition at line 63 of file DissipatingBody.h.