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 
79  std::vector<Dissipation::QuantityEntry> __orbit_entries;
80 
81  std::valarray<double>
87  __orbit_power,
88 
92 
100  std::vector<Eigen::Vector3d> __orbit_torque,
101 
108 
111 
114  std::valarray<Eigen::VectorXd> __above_lock_fractions;
115 
118  double normalize_torques(
120  double companion_mass,
121 
123  double semimajor,
124 
126  double orbital_frequency
127  );
128 
131  void collect_orbit_rates(
133  double orbital_frequency,
134 
137  double torque_norm
138  );
139 
142  void calculate_orbit_rate_corrections();
143 
149  void angular_momentum_transfer(
151  const DissipatingZone &outer_zone,
152 
154  const DissipatingZone &inner_zone,
155 
158  Eigen::Vector3d &outer_angmom_gain,
159 
162  Eigen::Vector3d &inner_angmom_gain,
163 
167  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
168 
172  bool with_respect_to_outer=false
173  ) const;
174 
182  Eigen::Vector3d angular_momentum_transfer_from_top(
185  unsigned zone_index,
186 
188  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
189 
193  bool with_respect_to_outer=false
194  ) const;
195 
202  Eigen::Vector3d angular_momentum_transfer_from_bottom(
205  unsigned zone_index,
206 
208  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
209 
213  bool with_respect_to_inner=false
214  ) const;
215 
222  Eigen::Vector3d angular_momentum_transfer_to_zone(
224  unsigned zone_index,
225 
227  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
228 
230  int deriv_zone=0
231  ) const;
232 
234  void calculate_nontidal_torques();
235 
237  void correct_orbit_power(
239  Eigen::VectorXd &above_lock_fractions_age_deriv,
240 
243  Eigen::VectorXd &above_lock_fractions_semimajor_deriv,
244 
247  Eigen::VectorXd &above_lock_fractions_eccentricity_deriv,
248 
251  Eigen::VectorXd &above_lock_fractions_radius_deriv
252  );
253 
255  void correct_orbit_torque(
257  std::valarray<Eigen::VectorXd> &above_lock_fractions
258  );
259  public:
261  DissipatingBody();
262 
267  virtual void configure(
269  bool initialize,
270 
272  double age,
273 
275  double companion_mass,
276 
278  double semimajor,
279 
281  double eccentricity,
282 
285  const double *spin_angmom,
286 
289  const double *inclination = NULL,
290 
294  const double *periapsis = NULL,
295 
298  bool locked_surface = false,
299 
303  bool zero_outer_inclination = false,
304 
308  bool zero_outer_periapsis = false
309  );
310 
313  void lock_zone_spin(unsigned zone_index,
314  int orbital_frequency_multiplier,
315  int spin_frequency_multiplier)
316  {
317  zone(zone_index).set_lock(orbital_frequency_multiplier,
318  spin_frequency_multiplier);
319  ++__num_locked_zones;
320  }
321 
325  unsigned zone_index,
326 
329  short direction
330  )
331  {
332  if(direction == 0) zone(zone_index).release_lock();
333  else zone(zone_index).release_lock(direction);
334  --__num_locked_zones;
335  }
336 
338  unsigned number_locked_zones() const {return __num_locked_zones;}
339 
350  Eigen::Vector3d nontidal_torque(
352  unsigned zone_index,
353 
355  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
356 
368  int deriv_zone=0
369  ) const;
370 
375  const Eigen::Vector3d &tidal_torque(
377  unsigned zone_index,
378 
382  bool above,
383 
385  Dissipation::QuantityEntry entry=Dissipation::NO_DERIV
386  ) const
387  {
388  assert(zone_index<number_zones());
389 
390  return (above
391  ? __tidal_torques_above
392  : __tidal_torques_below)[zone_index][entry];
393  }
394 
396  double tidal_power(
398  unsigned zone_index,
399 
403  bool above,
404 
406  Dissipation::QuantityEntry entry=Dissipation::NO_DERIV
407  ) const;
408 
411  void set_above_lock_fractions(
415  std::valarray<Eigen::VectorXd> &above_lock_fractions
416  );
417 
424  double tidal_orbit_power(
426  Dissipation::QuantityEntry entry=Dissipation::NO_DERIV,
427 
430  unsigned deriv_zone_index=0,
431 
434  const Eigen::VectorXd &
435  above_lock_fraction_deriv=Eigen::VectorXd()
436  ) const;
437 
443  Eigen::Vector3d tidal_orbit_torque(
445  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
446 
449  unsigned deriv_zone_index=0,
450 
453  const Eigen::VectorXd &
454  above_lock_fraction_deriv=Eigen::VectorXd()
455  ) const;
456 
460  Eigen::Vector3d tidal_orbit_torque(
462  const DissipatingZone &reference_zone,
463 
465  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
466 
469  unsigned deriv_zone_index=0,
470 
473  const Eigen::VectorXd &above_lock_fraction_deriv=
474  Eigen::VectorXd()
475  ) const;
476 
478  virtual unsigned number_zones() const =0;
479 
481  virtual const DissipatingZone &zone(
485  unsigned zone_index
486  ) const=0;
487 
489  virtual DissipatingZone &zone(
493  unsigned zone_index
494  )=0;
495 
501  virtual Eigen::Vector3d angular_momentum_coupling(
504  unsigned top_zone_index,
505 
507  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV,
508 
512  bool with_respect_to_top=false
513  ) const =0;
514 
520  virtual double angular_momentum_loss(
522  Dissipation::QuantityEntry deriv=Dissipation::NO_DERIV
523  ) const =0;
524 
526  double radius(
528  int deriv_order=0
529  ) const
530  {return zone(0).outer_radius(deriv_order);}
531 
533  double mass() const
534  {return zone(0).outer_mass(Dissipation::NO_DERIV);}
535 
537  double spin_frequency() const {return zone(0).spin_frequency();}
538 
544  double surface_lock_frequency() const
545  {return __surface_lock_frequency;}
546 
548  void set_surface_lock_frequency(double frequency)
549  {__surface_lock_frequency=frequency;}
550 
552  virtual void add_to_evolution();
553 
555  virtual void rewind_evolution(
557  unsigned nsteps
558  );
559 
561  virtual void reset_evolution();
562 
567  virtual CombinedStoppingCondition *stopping_conditions(
569  BinarySystem &system,
570 
572  bool primary
573  );
574 
576  virtual void spin_jumped() {zone(0).spin_jumped();}
577 
581  virtual void reached_critical_age(double) {assert(false);}
582 
585  virtual double next_stop_age() const {return Core::Inf;}
586 
588  virtual void change_expansion_order(
590  unsigned new_expansion_order,
591 
593  BinarySystem &system,
594 
596  bool primary
597  );
598 
600  virtual ~DissipatingBody() {}
601  }; //End DissipatingBody class.
602 
603 } //End Evolve namespace.
604 #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.
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:57
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.