9 #ifndef __DISSIPATING_BODY_H 10 #define __DISSIPATING_BODY_H 12 #include "../Core/SharedLibraryExportMacros.h" 15 #include "../Core/OrbitalExpressions.h" 16 #include "../Core/AstronomicalConstants.h" 17 #include "../Core/Common.h" 49 __dorbital_frequency_da,
54 std::valarray< std::valarray<Eigen::Vector3d> >
63 __angular_momentum_transfer,
69 __tidal_torques_above,
100 std::vector<Eigen::Vector3d> __orbit_torque,
118 double normalize_torques(
120 double companion_mass,
126 double orbital_frequency
131 void collect_orbit_rates(
133 double orbital_frequency,
142 void calculate_orbit_rate_corrections();
149 void angular_momentum_transfer(
158 Eigen::Vector3d &outer_angmom_gain,
162 Eigen::Vector3d &inner_angmom_gain,
172 bool with_respect_to_outer=
false 182 Eigen::Vector3d angular_momentum_transfer_from_top(
193 bool with_respect_to_outer=
false 202 Eigen::Vector3d angular_momentum_transfer_from_bottom(
213 bool with_respect_to_inner=
false 222 Eigen::Vector3d angular_momentum_transfer_to_zone(
234 void calculate_nontidal_torques();
237 void correct_orbit_power(
239 Eigen::VectorXd &above_lock_fractions_age_deriv,
243 Eigen::VectorXd &above_lock_fractions_semimajor_deriv,
247 Eigen::VectorXd &above_lock_fractions_eccentricity_deriv,
251 Eigen::VectorXd &above_lock_fractions_radius_deriv
255 void correct_orbit_torque(
257 std::valarray<Eigen::VectorXd> &above_lock_fractions
267 virtual void configure(
275 double companion_mass,
285 const double *spin_angmom,
289 const double *inclination = NULL,
294 const double *periapsis = NULL,
298 bool locked_surface =
false,
303 bool zero_outer_inclination =
false,
308 bool zero_outer_periapsis =
false 314 int orbital_frequency_multiplier,
315 int spin_frequency_multiplier)
317 zone(zone_index).set_lock(orbital_frequency_multiplier,
318 spin_frequency_multiplier);
319 ++__num_locked_zones;
332 if(direction == 0) zone(zone_index).release_lock();
333 else zone(zone_index).release_lock(direction);
334 --__num_locked_zones;
350 Eigen::Vector3d nontidal_torque(
388 assert(zone_index<number_zones());
391 ? __tidal_torques_above
392 : __tidal_torques_below)[zone_index][entry];
411 void set_above_lock_fractions(
415 std::valarray<Eigen::VectorXd> &above_lock_fractions
424 double tidal_orbit_power(
430 unsigned deriv_zone_index=0,
434 const Eigen::VectorXd &
435 above_lock_fraction_deriv=Eigen::VectorXd()
443 Eigen::Vector3d tidal_orbit_torque(
449 unsigned deriv_zone_index=0,
453 const Eigen::VectorXd &
454 above_lock_fraction_deriv=Eigen::VectorXd()
460 Eigen::Vector3d tidal_orbit_torque(
469 unsigned deriv_zone_index=0,
473 const Eigen::VectorXd &above_lock_fraction_deriv=
478 virtual unsigned number_zones()
const =0;
501 virtual Eigen::Vector3d angular_momentum_coupling(
504 unsigned top_zone_index,
512 bool with_respect_to_top=
false 520 virtual double angular_momentum_loss(
530 {
return zone(0).outer_radius(deriv_order);}
545 {
return __surface_lock_frequency;}
549 {__surface_lock_frequency=frequency;}
552 virtual void add_to_evolution();
555 virtual void rewind_evolution(
561 virtual void reset_evolution();
588 virtual void change_expansion_order(
590 unsigned new_expansion_order,
unsigned number_locked_zones() const
The number of zones currently in a spin-orbit lock.
virtual ~DissipatingBody()
Virtual destructor.
Declares a class representing one zone of a body dissipative to tidal distortions.
A base class for any body contributing to tidal dissipation.
void unlock_zone_spin(unsigned zone_index, short direction)
Releases the given zone from a spin-orbit lock.
virtual void reached_critical_age(double)
Change the body as necessary at the given age.
double surface_lock_frequency() const
Angular velocity of the surface zone when locked (assumed constant).
Orientations of zones of bodies in a binary system.
double mass() const
The mass of the body (constant with age).
double radius(int deriv_order=0) const
The current radius or its derivative with age of the body.
virtual void spin_jumped()
Notifies the body that its spin just discontinously jumped.
A layer of a system body for which the tidal bulge is not exactly in phase with the tidal potential...
virtual double next_stop_age() const
The next age when the evolution needs to be stopped for a change in one of the bodies.
std::valarray< std::valarray< Eigen::Vector3d > > __tidal_torques_below
Tidal torques and their derivatives on each zone for spin frequency approaching potential lock from b...
std::vector< Eigen::Vector3d > __orbit_torque_correction
Corrections to __orbit_torque_below (undifferentiated) if single zones switch to above.
NO_DERIV
The quantity itself, undifferentiated.
std::vector< Dissipation::QuantityEntry > __orbit_entries
The quantities w.r.t. which derivatives of the orbit energy gain and torque are pre-calculated.
void lock_zone_spin(unsigned zone_index, int orbital_frequency_multiplier, int spin_frequency_multiplier)
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()).
void set_surface_lock_frequency(double frequency)
Sets the frequency at which the surface is locked (if any).
double __surface_lock_frequency
The frequency at which the surface is locked (if any).
A class combining the the outputs of multiple stopping conditions.
std::valarray< Eigen::VectorXd > __above_lock_fractions
The fractional contribution of the above the lock rates for locked zones and their derivatives...
unsigned __num_locked_zones
The number of zones currently in a spin-orbit lock.
Describes a system of two bodies orbiting each other.
std::valarray< double > __orbit_power_correction
Corrections to __orbit_power_below (undifferentiated) if single zones switch to above.
double spin_frequency() const
The surface spin freuqency of the body.