Planetary Orbital Evolution due to Tides
Orbital evolution of two objects experiencing tides
DissipatingBody.h
Go to the documentation of this file.
1 
9 #ifndef __DISSIPATING_BODY_H
10 #define __DISSIPATING_BODY_H
11 
12 #include "../Core/SharedLibraryExportMacros.h"
13 
14 #include "DissipatingZone.h"
15 #include "../Core/OrbitalExpressions.h"
16 #include "../Core/AstronomicalConstants.h"
17 #include "../Core/Common.h"
18 #include <valarray>
19 #include <cassert>
20 
21 namespace Evolve {
22 
23  class BinarySystem;
24 
39  class LIB_PUBLIC DissipatingBody {
40  private:
41  double
43  __power_norm,
44 
46  __orbital_frequency,
47 
49  __dorbital_frequency_da,
50 
53 
54  std::valarray< std::valarray<Eigen::Vector3d> >
63  __angular_momentum_transfer,
64 
69  __tidal_torques_above,
70 
76 
80  std::vector<Dissipation::QuantityEntry> __orbit_entries;
81 
82  std::valarray<double>
88  __orbit_power,
89 
93 
101  std::vector<Eigen::Vector3d> __orbit_torque,
102 
109 
112 
115  std::valarray<Eigen::VectorXd> __above_lock_fractions;
116 
119  double normalize_torques(
121  double companion_mass,
122 
124  double semimajor,
125 
127  double orbital_frequency
128  );
129 
132  void collect_orbit_rates(
134  double orbital_frequency,
135 
138  double torque_norm
139  );
140 
143  void calculate_orbit_rate_corrections();
144 
150  void angular_momentum_transfer(
152  const DissipatingZone &outer_zone,
153 
155  const DissipatingZone &inner_zone,
156 
159  Eigen::Vector3d &outer_angmom_gain,
160 
163  Eigen::Vector3d &inner_angmom_gain,
164 
168  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
169 
173  bool with_respect_to_outer=false
174  ) const;
175 
183  Eigen::Vector3d angular_momentum_transfer_from_top(
186  unsigned zone_index,
187 
189  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
190 
194  bool with_respect_to_outer=false
195  ) const;
196 
203  Eigen::Vector3d angular_momentum_transfer_from_bottom(
206  unsigned zone_index,
207 
209  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
210 
214  bool with_respect_to_inner=false
215  ) const;
216 
223  Eigen::Vector3d angular_momentum_transfer_to_zone(
225  unsigned zone_index,
226 
228  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
229 
231  int deriv_zone=0
232  ) const;
233 
235  void calculate_nontidal_torques();
236 
238  void correct_orbit_power(
240  Eigen::VectorXd &above_lock_fractions_age_deriv,
241 
244  Eigen::VectorXd &above_lock_fractions_semimajor_deriv,
245 
248  Eigen::VectorXd &above_lock_fractions_eccentricity_deriv,
249 
252  Eigen::VectorXd &above_lock_fractions_radius_deriv
253  );
254 
256  void correct_orbit_torque(
258  std::valarray<Eigen::VectorXd> &above_lock_fractions
259  );
260  public:
262  DissipatingBody();
263 
268  virtual void configure(
270  bool initialize,
271 
273  double age,
274 
276  double companion_mass,
277 
279  double semimajor,
280 
282  double eccentricity,
283 
286  const double *spin_angmom,
287 
290  const double *inclination = NULL,
291 
295  const double *periapsis = NULL,
296 
299  bool locked_surface = false,
300 
304  bool zero_outer_inclination = false,
305 
309  bool zero_outer_periapsis = false
310  );
311 
314  void lock_zone_spin(unsigned zone_index,
315  int orbital_frequency_multiplier,
316  int spin_frequency_multiplier)
317  {
318  zone(zone_index).set_lock(orbital_frequency_multiplier,
319  spin_frequency_multiplier);
320  ++__num_locked_zones;
321  }
322 
326  unsigned zone_index,
327 
330  short direction
331  )
332  {
333  if(direction == 0) zone(zone_index).release_lock();
334  else zone(zone_index).release_lock(direction);
335  --__num_locked_zones;
336  }
337 
339  unsigned number_locked_zones() const {return __num_locked_zones;}
340 
351  Eigen::Vector3d nontidal_torque(
353  unsigned zone_index,
354 
356  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
357 
369  int deriv_zone=0
370  ) const;
371 
376  const Eigen::Vector3d &tidal_torque(
378  unsigned zone_index,
379 
383  bool above,
384 
386  Dissipation::QuantityEntry entry=Dissipation::NO_DERIV
387  ) const
388  {
389  assert(zone_index<number_zones());
390 
391  return (above
392  ? __tidal_torques_above
393  : __tidal_torques_below)[zone_index][entry];
394  }
395 
397  double tidal_power(
399  unsigned zone_index,
400 
404  bool above,
405 
407  Dissipation::QuantityEntry entry=Dissipation::NO_DERIV
408  ) const;
409 
412  void set_above_lock_fractions(
416  std::valarray<Eigen::VectorXd> &above_lock_fractions
417  );
418 
425  double tidal_orbit_power(
427  Dissipation::QuantityEntry entry=Dissipation::NO_DERIV,
428 
431  unsigned deriv_zone_index=0,
432 
435  const Eigen::VectorXd &
436  above_lock_fraction_deriv=Eigen::VectorXd()
437  ) const;
438 
444  Eigen::Vector3d tidal_orbit_torque(
446  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
447 
450  unsigned deriv_zone_index=0,
451 
454  const Eigen::VectorXd &
455  above_lock_fraction_deriv=Eigen::VectorXd()
456  ) const;
457 
461  Eigen::Vector3d tidal_orbit_torque(
463  const DissipatingZone &reference_zone,
464 
466  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
467 
470  unsigned deriv_zone_index=0,
471 
474  const Eigen::VectorXd &above_lock_fraction_deriv=
475  Eigen::VectorXd()
476  ) const;
477 
479  virtual unsigned number_zones() const =0;
480 
482  virtual const DissipatingZone &zone(
486  unsigned zone_index
487  ) const=0;
488 
490  virtual DissipatingZone &zone(
494  unsigned zone_index
495  )=0;
496 
502  virtual Eigen::Vector3d angular_momentum_coupling(
505  unsigned top_zone_index,
506 
508  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
509 
513  bool with_respect_to_top=false
514  ) const =0;
515 
521  virtual double angular_momentum_loss(
523  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV
524  ) const =0;
525 
527  double radius(
529  int deriv_order=0
530  ) const
531  {return zone(0).outer_radius(deriv_order);}
532 
534  double mass() const
535  {return zone(0).outer_mass(Dissipation::NO_DERIV);}
536 
538  double spin_frequency() const {return zone(0).spin_frequency();}
539 
545  double surface_lock_frequency() const
546  {return __surface_lock_frequency;}
547 
549  void set_surface_lock_frequency(double frequency)
550  {__surface_lock_frequency=frequency;}
551 
553  virtual void add_to_evolution();
554 
556  virtual void rewind_evolution(
558  unsigned nsteps
559  );
560 
562  virtual void reset_evolution();
563 
568  virtual CombinedStoppingCondition *stopping_conditions(
570  BinarySystem &system,
571 
573  bool primary
574  );
575 
577  virtual void spin_jumped() {zone(0).spin_jumped();}
578 
582  virtual void reached_critical_age(double) {assert(false);}
583 
586  virtual double next_stop_age() const {return Core::Inf;}
587 
589  virtual void change_e_order(
591  unsigned new_e_order,
592 
594  BinarySystem &system,
595 
597  bool primary
598  );
599 
601  virtual ~DissipatingBody() {}
602  }; //End DissipatingBody class.
603 
604 } //End Evolve namespace.
605 #endif
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 or err...
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.
Definition: BinarySystem.h:56
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.