Planetary Orbital Evolution due to Tides
Orbital evolution of two objects experiencing tides
DissipatingBody.cpp
1 #define BUILDING_LIBRARY
2 #include "DissipatingBody.h"
3 
4 namespace Evolve {
5 
6  double DissipatingBody::normalize_torques(double companion_mass,
7  double semimajor,
8  double orbital_frequency)
9  {
10  double torque_norm = (
11  std::pow(companion_mass / std::pow(semimajor, 3), 2)
12  *
13  std::pow(radius(), 5)
14  *
16  *
18  *
20  *
22  /
24  );
25  double dangular_velocity_da = -1.5 * orbital_frequency / semimajor;
26  for(unsigned zone_index = 0; zone_index < number_zones(); ++zone_index) {
27  bool above = false;
28  const DissipatingZone &current_zone = zone(zone_index);
29  do {
30  std::valarray<Eigen::Vector3d> &tidal_torque = (
31  above
34  )[zone_index];
35  tidal_torque[Dissipation::ORBITAL_FREQUENCY] += (
36  4.0
37  /
38  orbital_frequency
39  *
40  tidal_torque[Dissipation::NO_DERIV]
41  );
42  tidal_torque[Dissipation::RADIUS] = (
43  5.0 / radius() * tidal_torque[Dissipation::NO_DERIV]
44  );
45  tidal_torque[Dissipation::MOMENT_OF_INERTIA] = (
46  -current_zone.spin_frequency()
47  /
48  current_zone.moment_of_inertia()
49  *
50  tidal_torque[Dissipation::SPIN_FREQUENCY]
51  );
52  tidal_torque[Dissipation::SPIN_ANGMOM] = (
53  tidal_torque[Dissipation::SPIN_FREQUENCY]
54  /
55  current_zone.moment_of_inertia()
56  );
57  tidal_torque[Dissipation::SEMIMAJOR] = (
58  dangular_velocity_da
59  *
60  tidal_torque[Dissipation::ORBITAL_FREQUENCY]
61  );
62  above = !above;
63  } while (above);
64  }
65 
66  assert(__tidal_torques_above.size() == __tidal_torques_below.size());
67  for(unsigned i = 0; i < __tidal_torques_above.size(); ++i) {
68  assert(__tidal_torques_above[i].size()
69  ==
70  __tidal_torques_below[i].size());
71  for(unsigned j = 0; j < __tidal_torques_above[i].size(); ++j) {
72  __tidal_torques_above[i][j] *= torque_norm;
73  __tidal_torques_below[i][j] *= torque_norm;
74  }
75  }
76  return torque_norm;
77  }
78 
79  void DissipatingBody::collect_orbit_rates(double orbital_frequency,
80  double torque_norm)
81  {
82  __power_norm = torque_norm * orbital_frequency;
83  DissipatingZone &surface_zone = zone(0);
84  unsigned invalid_ind = __orbit_entries.size(),
85  no_deriv_ind = invalid_ind,
86  orbital_freq_deriv_ind = invalid_ind,
87  radius_deriv_ind = invalid_ind,
88  semimajor_deriv_ind = invalid_ind;
89  for(unsigned i=0; i<__orbit_entries.size(); ++i)
90  switch(__orbit_entries[i]) {
91  case Dissipation::NO_DERIV : no_deriv_ind=i; break;
92  case Dissipation::ORBITAL_FREQUENCY : orbital_freq_deriv_ind=i;
93  break;
94  case Dissipation::RADIUS : radius_deriv_ind=i; break;
95  case Dissipation::SEMIMAJOR : semimajor_deriv_ind=i; break;
96  default:;
97  }
98  __orbit_power = 0;
99  for(
100  unsigned zone_index = 0;
101  zone_index < number_zones();
102  ++zone_index
103  ) {
104  DissipatingZone &current_zone = zone(zone_index);
105  std::valarray<Eigen::Vector3d> &tidal_torque = (
106  __tidal_torques_below[zone_index]
107  );
108  for(
109  unsigned deriv_ind = 0;
110  deriv_ind < __orbit_entries.size();
111  ++deriv_ind
112  ) {
113  Dissipation::QuantityEntry entry = __orbit_entries[deriv_ind];
115  __orbit_power[deriv_ind] -= (
116  current_zone.tidal_power(false, entry)
117  );
118  if(zone_index) {
120  current_zone,
121  surface_zone,
122  tidal_torque[entry]
123  );
124  } else __orbit_torque[entry] = -tidal_torque[entry];
125  }
126  if(zone_index) {
128  zone_to_zone_transform(current_zone,
129  surface_zone,
130  tidal_torque[Dissipation::NO_DERIV],
132  false)
133  );
135  zone_to_zone_transform(current_zone,
136  surface_zone,
137  tidal_torque[Dissipation::NO_DERIV],
139  false)
140  );
141  }
142  }
144  __orbit_power[orbital_freq_deriv_ind] += (
145  5.0 / orbital_frequency * __orbit_power[no_deriv_ind]
146  );
147  __orbit_power[radius_deriv_ind] += (
148  5.0 / radius() * __orbit_power[no_deriv_ind]
149  );
150  __orbit_power[semimajor_deriv_ind] = (
152  *
153  __orbit_power[orbital_freq_deriv_ind]
154  );
155  }
156 
158  {
161  unsigned correction_index = 0;
162  for(unsigned zone_index=0; zone_index<number_zones(); ++zone_index)
163  if(zone(zone_index).locked()) {
164  DissipatingZone &this_zone=zone(zone_index);
165  __orbit_power_correction[correction_index] = (
166  tidal_power(zone_index, false)
167  -
168  tidal_power(zone_index, true)
169  );
170  __orbit_torque_correction[correction_index] = (
171  zone_to_zone_transform(this_zone, zone(0),
172  tidal_torque(zone_index, false)
173  -
174  tidal_torque(zone_index, true))
175  );
176  ++correction_index;
177  }
178  }
179 
181  const DissipatingZone &outer_zone,
182  const DissipatingZone &inner_zone,
183  Eigen::Vector3d &outer_angmom_gain,
184  Eigen::Vector3d &inner_angmom_gain,
185  Dissipation::QuantityEntry deriv,
186  bool with_respect_to_outer
187  ) const
188  {
189  assert(deriv == Dissipation::NO_DERIV
190  ||
191  deriv == Dissipation::INCLINATION
192  ||
193  deriv == Dissipation::PERIAPSIS);
194 
195  double dm_dt = inner_zone.outer_mass(1),
196  lost_spin = (dm_dt >= 0
197  ? outer_zone
198  : inner_zone).spin_frequency(),
199  angmom_transfer = (
200  -2.0 / 3.0
201  *
202  std::pow(inner_zone.outer_radius(), 2)
203  *
204  lost_spin * dm_dt
205  );
206  if(dm_dt > 0) {
207  outer_angmom_gain[0] = outer_angmom_gain[1] = 0;
208  outer_angmom_gain[2] = (deriv == Dissipation::NO_DERIV
209  ? angmom_transfer
210  : 0);
211  inner_angmom_gain = -zone_to_zone_transform(
212  outer_zone,
213  inner_zone,
214  Eigen::Vector3d(0, 0, angmom_transfer),
215  deriv,
216  with_respect_to_outer
217  );
218  } else {
219  inner_angmom_gain[0] = inner_angmom_gain[1] = 0;
220  inner_angmom_gain[2] = (deriv == Dissipation::NO_DERIV
221  ? -angmom_transfer
222  : 0);
223  outer_angmom_gain=zone_to_zone_transform(
224  inner_zone,
225  outer_zone,
226  Eigen::Vector3d(0, 0, angmom_transfer),
227  deriv,
228  !with_respect_to_outer
229  );
230  }
231  }
232 
234  unsigned zone_index,
235  Dissipation::QuantityEntry deriv,
236  bool with_respect_to_outer
237  ) const
238  {
239  assert(zone_index > 0);
240 
241  const DissipatingZone &this_zone = zone(zone_index),
242  &zone_above = zone(zone_index - 1);
243  double scaling = Core::NaN,
244  dm_dt = this_zone.outer_mass(1);
245  if(deriv == Dissipation::NO_DERIV)
246  scaling = 1;
247  else if(deriv == Dissipation::AGE) {
248  scaling = (
249  2.0 * this_zone.outer_radius(1) / this_zone.outer_radius(0)
250  +
251  this_zone.outer_mass(2) / dm_dt
252  );
253  } else if(
255  ||
257  ||
258  deriv == Dissipation::SPIN_ANGMOM
259  ) {
260  if(
261  (with_respect_to_outer && dm_dt <= 0)
262  ||
263  (!with_respect_to_outer && dm_dt>=0)
264  )
265  return Eigen::Vector3d(0, 0, 0);
266  else if(deriv == Dissipation::SPIN_FREQUENCY)
267  scaling = 1.0 / (with_respect_to_outer
268  ? zone_above.spin_frequency()
269  : this_zone.spin_frequency());
270  else if(deriv == Dissipation::MOMENT_OF_INERTIA)
271  scaling = -1.0 / (with_respect_to_outer
272  ? zone_above.moment_of_inertia()
273  : this_zone.moment_of_inertia());
274  else scaling = 1.0 / (with_respect_to_outer
275  ? zone_above.angular_momentum()
276  : this_zone.angular_momentum());
277  } else if(
278  deriv == Dissipation::INCLINATION
279  ||
280  deriv == Dissipation::PERIAPSIS
281  ) {
282  if(dm_dt <= 0)
283  return Eigen::Vector3d(0, 0, 0);
284  else {
285  Eigen::Vector3d dummy, result;
286  angular_momentum_transfer(zone_above,
287  this_zone,
288  dummy,
289  result,
290  deriv,
291  with_respect_to_outer);
292  return result;
293  }
294  } else
295  assert(false);
296  return scaling * __angular_momentum_transfer[zone_index - 1][1];
297  }
298 
300  unsigned zone_index,
301  Dissipation::QuantityEntry deriv,
302  bool with_respect_to_inner
303  ) const
304  {
305  assert(zone_index < number_zones() - 1);
306 
307  const DissipatingZone &this_zone = zone(zone_index),
308  &zone_below = zone(zone_index + 1);
309  double scaling = Core::NaN,
310  dm_dt = zone_below.outer_mass(1);
311  if(deriv == Dissipation::NO_DERIV)
312  scaling = 1;
313  else if(deriv == Dissipation::AGE) {
314  scaling = (
315  2.0 * zone_below.outer_radius(1) / zone_below.outer_radius(0)
316  +
317  zone_below.outer_mass(2) / dm_dt
318  );
319  } else if(
321  ||
323  ||
324  deriv == Dissipation::SPIN_ANGMOM
325  ) {
326  if(
327  (with_respect_to_inner && dm_dt >= 0)
328  ||
329  (!with_respect_to_inner && dm_dt <= 0)
330  )
331  return Eigen::Vector3d(0, 0, 0);
332  else if(deriv == Dissipation::SPIN_FREQUENCY)
333  scaling = 1.0 / (with_respect_to_inner
334  ? zone_below.spin_frequency()
335  : this_zone.spin_frequency());
336  else if(deriv == Dissipation::MOMENT_OF_INERTIA)
337  scaling = -1.0 / (with_respect_to_inner
338  ? zone_below.moment_of_inertia()
339  : this_zone.moment_of_inertia());
340  else scaling = 1.0 / (with_respect_to_inner
341  ? zone_below.angular_momentum()
342  : this_zone.angular_momentum());
343  } else if(
344  deriv == Dissipation::INCLINATION
345  ||
346  deriv == Dissipation::PERIAPSIS
347  ) {
348  if(dm_dt >= 0)
349  return Eigen::Vector3d(0, 0, 0);
350  else {
351  Eigen::Vector3d dummy, result;
352  angular_momentum_transfer(this_zone,
353  zone_below,
354  result,
355  dummy,
356  deriv,
357  !with_respect_to_inner);
358  }
359  } else
360  assert(false);
361 
362  return scaling * __angular_momentum_transfer[zone_index][0];
363  }
364 
366  unsigned zone_index,
367  Dissipation::QuantityEntry deriv,
368  int deriv_zone
369  ) const
370  {
371 
372  assert(deriv < Dissipation::NUM_DERIVATIVES);
373  if(
375  ||
377  ||
378  deriv == Dissipation::RADIUS
379  ||
380  deriv == Dissipation::SEMIMAJOR
381  )
382  return Eigen::Vector3d(0,0,0);
383  Eigen::Vector3d result(0, 0, 0);
384  if(
385  zone_index > 0
386  &&
387  (
388  deriv == Dissipation::NO_DERIV
389  ||
390  deriv == Dissipation::AGE
391  ||
392  deriv_zone <= 0
393  )
394  )
395  result = angular_momentum_transfer_from_top(zone_index,
396  deriv,
397  deriv_zone < 0);
398  if(
399  zone_index < number_zones() - 1
400  &&
401  (
402  deriv == Dissipation::NO_DERIV
403  ||
404  deriv == Dissipation::AGE
405  ||
406  deriv_zone >= 0
407  )
408  )
409  result += angular_momentum_transfer_from_bottom(zone_index,
410  deriv,
411  deriv_zone > 0);
412  return result;
413  }
414 
416  __orbit_entries(6),
418  __orbit_torque(Dissipation::NUM_ENTRIES),
421  {
428  }
429 
431  Eigen::VectorXd &above_lock_fractions_age_deriv,
432  Eigen::VectorXd &above_lock_fractions_semimajor_deriv,
433  Eigen::VectorXd &above_lock_fractions_eccentricity_deriv,
434  Eigen::VectorXd &above_lock_fractions_radius_deriv
435  )
436  {
437  unsigned correction_index = 0;
438  for(
439  unsigned zone_index = 0;
440  zone_index < number_zones();
441  ++zone_index
442  ) {
443  if(zone(zone_index).locked()) {
444  unsigned locked_zone_index = (
445  zone(zone_index).locked_zone_index()
446  );
447  __orbit_power_correction[correction_index] = (
448  tidal_power(zone_index, false)
449  -
450  tidal_power(zone_index, true)
451  );
452  for(
453  unsigned entry_ind = 0;
454  entry_ind < __orbit_entries.size();
455  ++entry_ind
456  ) {
457  Dissipation::QuantityEntry entry = (
458  __orbit_entries[entry_ind]
459  );
460  __orbit_power[entry_ind] += (
463  [locked_zone_index]
464  *
465  (
466  tidal_power(zone_index, false, entry)
467  -
468  tidal_power(zone_index, true, entry)
469  )
470  );
471  if(entry != Dissipation::NO_DERIV) {
472  double frac_deriv = Core::NaN;
473  if(entry == Dissipation::AGE)
474  frac_deriv = (above_lock_fractions_age_deriv
475  [locked_zone_index]);
476  else if(entry == Dissipation::ORBITAL_FREQUENCY)
477  frac_deriv = (
478  above_lock_fractions_semimajor_deriv
479  [locked_zone_index]
480  /
482  );
483  else if(entry == Dissipation::ECCENTRICITY)
484  frac_deriv = (
485  above_lock_fractions_eccentricity_deriv
486  [locked_zone_index]
487  /
489  );
490  else if(entry == Dissipation::RADIUS)
491  frac_deriv = (above_lock_fractions_radius_deriv
492  [locked_zone_index]);
493  else if(entry == Dissipation::SEMIMAJOR)
494  frac_deriv = (above_lock_fractions_semimajor_deriv
495  [locked_zone_index]);
496  else
497  assert(false);
498 
499  __orbit_power[entry_ind] += (
500  frac_deriv
501  *
502  __orbit_power_correction[correction_index]
503  );
504  }
505  }
506  ++correction_index;
507  }
508  }
509  }
510 
512  std::valarray<Eigen::VectorXd> &above_lock_fractions
513  )
514  {
515  DissipatingZone &surface_zone = zone(0);
516  unsigned correction_index = 0;
517  for(
518  unsigned zone_index = 0;
519  zone_index < number_zones();
520  ++zone_index
521  ) {
522  if(zone(zone_index).locked()) {
523  unsigned locked_zone_index = (
524  zone(zone_index).locked_zone_index()
525  );
526  DissipatingZone &this_zone = zone(zone_index);
527  for(
528  int entry_int = Dissipation::NO_DERIV;
529  entry_int < Dissipation::NUM_ENTRIES;
530  ++entry_int
531  ) {
532  Dissipation::QuantityEntry entry = (
533  static_cast<Dissipation::QuantityEntry>(entry_int)
534  );
535  Eigen::Vector3d correction = (
536  above_lock_fractions
538  [locked_zone_index]
539  *
540  (
541  tidal_torque(zone_index, false, entry)
542  -
543  tidal_torque(zone_index, true, entry)
544  )
545  );
546  if(zone_index) {
548  this_zone,
549  surface_zone,
550  correction
551  );
552  if(
553  entry == Dissipation::INCLINATION
554  ||
555  entry == Dissipation::PERIAPSIS
556  )
558  this_zone,
559  surface_zone,
560  (
561  above_lock_fractions
563  [locked_zone_index]
564  *
566  [correction_index]
567  ),
568  static_cast<Dissipation::QuantityEntry>(
569  entry
570  )
571  );
572  } else
573  __orbit_torque[entry] += correction;
574 
575  if(entry != Dissipation::NO_DERIV)
576  __orbit_torque[entry] += (
577  above_lock_fractions[entry][locked_zone_index]
578  *
579  __orbit_torque_correction[correction_index]
580  );
581  }
582  ++correction_index;
583  }
584  }
585  }
586 
587  void DissipatingBody::configure(bool initialize,
588  double age,
589  double companion_mass,
590  double semimajor,
591  double eccentricity,
592  const double *spin_angmom,
593  const double *inclination,
594  const double *periapsis,
595  bool locked_surface,
596  bool zero_outer_inclination,
597  bool zero_outer_periapsis)
598  {
599  double orbital_angmom=Core::orbital_angular_momentum(
600  companion_mass,
601  mass(),
602  semimajor,
603  eccentricity
604  );
605  __orbital_frequency=Core::orbital_angular_velocity(
606  companion_mass,
607  mass(),
608  semimajor,
609  false
610  );
611  __dorbital_frequency_da=Core::orbital_angular_velocity(
612  companion_mass,
613  mass(),
614  semimajor,
615  true
616  );
617 
618  if(initialize) {
619  __num_locked_zones = 0;
620 #ifndef NDEBUG
621  std::cerr << "Initializing DissipatingBody" << std::endl;
622 #endif
623  }
624 
628  unsigned angmom_offset = (locked_surface ? 1 : 0);
629  for(
630  unsigned zone_index = 0;
631  zone_index < number_zones();
632  ++zone_index
633  ) {
634  DissipatingZone &current_zone = zone(zone_index);
635  double zone_inclination, zone_periapsis, zone_spin;
636  if(!inclination)
637  zone_inclination = 0;
638  else if(zero_outer_inclination)
639  zone_inclination = (zone_index
640  ? inclination[zone_index - 1]
641  : 0);
642  else
643  zone_inclination = inclination[zone_index];
644  if(!periapsis)
645  zone_periapsis = 0;
646  else if(zero_outer_periapsis)
647  zone_periapsis = (zone_index ? periapsis[zone_index-1] : 0);
648  else
649  zone_periapsis = periapsis[zone_index];
650  if(locked_surface && zone_index == 0)
651  zone_spin = surface_lock_frequency();
652  else if(current_zone.locked() && !initialize) {
653  zone_spin = Core::NaN;
654  ++angmom_offset;
655  } else
656  zone_spin = spin_angmom[zone_index - angmom_offset];
657 #ifndef NDEBUG
658  std::cerr << "Configuring zone " << zone_index << std::endl;
659 #endif
660  current_zone.configure(initialize,
661  age,
663  eccentricity,
664  orbital_angmom,
665  zone_spin,
666  zone_inclination,
667  zone_periapsis,
668  locked_surface && zone_index == 0);
669  }
670  for(
671  unsigned zone_index = 0;
672  zone_index < number_zones();
673  ++zone_index
674  ) {
675  DissipatingZone &current_zone = zone(zone_index);
676  if(zone_index < number_zones() - 1) {
677  __angular_momentum_transfer[zone_index].resize(2);
679  current_zone,
680  zone(zone_index + 1),
681  __angular_momentum_transfer[zone_index][0],
682  __angular_momentum_transfer[zone_index][1]
683  );
684  }
685  bool above = false;
686  do {
687  std::valarray<Eigen::Vector3d> &tidal_torque =
689  [zone_index];
690  tidal_torque.resize(Dissipation::NUM_ENTRIES);
691  for(
692  int torque_ind = Dissipation::NO_DERIV;
694  ++torque_ind
695  ) {
696  Dissipation::QuantityEntry entry = (
697  static_cast<Dissipation::QuantityEntry>(torque_ind)
698  );
699  tidal_torque[entry][0] = current_zone.tidal_torque_x(
700  above,
701  entry
702  );
703  tidal_torque[entry][1] = current_zone.tidal_torque_y(
704  above,
705  entry
706  );
707  tidal_torque[entry][2] = current_zone.tidal_torque_z(
708  above,
709  entry
710  );
711  //TODO: revive eccentricity derivative check
712  assert(
714  ||
715  !std::isnan(tidal_torque[entry].sum())
716  );
717  }
718  above = !above;
719  } while(above);
720  }
722  normalize_torques(companion_mass,
723  semimajor,
726  __above_lock_fractions.resize(0);
727  }
728 
730  unsigned zone_index,
731  Dissipation::QuantityEntry deriv,
732  int deriv_zone
733  ) const
734  {
735  assert(deriv < Dissipation::NUM_DERIVATIVES);
736  assert(zone_index < number_zones());
737  assert(static_cast<int>(zone_index) + deriv_zone >= 0);
738  assert(static_cast<int>(zone_index) + deriv_zone
739  <
740  static_cast<int>(number_zones()));
741  const DissipatingZone &this_zone = zone(zone_index);
742  Eigen::Vector3d result(0, 0, 0);
743  if(zone_index == 0 && deriv_zone == 0)
744  result[2] = -angular_momentum_loss(deriv);
745  result += angular_momentum_transfer_to_zone(zone_index,
746  deriv,
747  deriv_zone);
748  if(
749  zone_index < number_zones() - 1
750  &&
751  (!zone_specific(deriv) || deriv_zone >= 0)
752  )
753  result += angular_momentum_coupling(zone_index,
754  deriv,
755  deriv_zone == 0);
756  if(
757  zone_index > 0
758  &&
759  (!zone_specific(deriv) || deriv_zone<=0)
760  ) {
761  result -= zone_to_zone_transform(
762  zone(zone_index-1),
763  this_zone,
764  angular_momentum_coupling(zone_index - 1,
765  deriv,
766  deriv_zone < 0)
767  );
768  if(
769  deriv == Dissipation::INCLINATION
770  ||
771  deriv == Dissipation::PERIAPSIS
772  )
773  result -= zone_to_zone_transform(
774  zone(zone_index - 1),
775  this_zone,
776  angular_momentum_coupling(zone_index - 1),
777  deriv,
778  deriv_zone<0
779  );
780  }
781  return result;
782  }
783 
784  double DissipatingBody::tidal_power(unsigned zone_index,
785  bool above,
786  Dissipation::QuantityEntry entry) const
787  {
788  assert(zone_index < number_zones());
789 
790  const DissipatingZone &this_zone = zone(zone_index);
791  double result = (entry < Dissipation::END_DIMENSIONLESS_DERIV
792  ? __power_norm * this_zone.tidal_power(above, entry)
793  : 0.0);
794  if(
796  ||
797  entry == Dissipation::SEMIMAJOR
798  ) {
799  result += (5.0
800  /
802  *
804  *
805  this_zone.tidal_power(above));
806  if(entry == Dissipation::SEMIMAJOR)
807  result *= __dorbital_frequency_da;
808  }
809  else if(entry == Dissipation::RADIUS)
810  result += 5.0 / radius() * this_zone.tidal_power(above);
811  return result;
812  }
813 
815  std::valarray<Eigen::VectorXd> &above_lock_fractions
816  )
817  {
818  assert(above_lock_fractions.size() == Dissipation::NUM_DERIVATIVES);
819 
821  for(unsigned i = 0; i < Dissipation::NUM_DERIVATIVES; ++i) {
822  assert(above_lock_fractions[i].size() >= __num_locked_zones);
823 
824  __above_lock_fractions[i].resize(above_lock_fractions[i].size());
825  }
826  __above_lock_fractions = above_lock_fractions;
828  above_lock_fractions[Dissipation::AGE],
829  above_lock_fractions[Dissipation::SEMIMAJOR],
830  above_lock_fractions[Dissipation::ECCENTRICITY],
831  above_lock_fractions[Dissipation::RADIUS]
832  );
833  correct_orbit_torque(above_lock_fractions);
834  }
835 
837  Dissipation::QuantityEntry entry,
838  unsigned deriv_zone_index,
839  const Eigen::VectorXd &above_lock_fraction_deriv
840  ) const
841  {
842  const DissipatingZone &deriv_zone = zone(deriv_zone_index);
843  double result = Core::NaN;
844  if(
845  !zone_specific(entry)
846  ||
847  !deriv_zone.locked()
848  ||
849  __above_lock_fractions.size() == 0
850  ) {
851  if(zone_specific(entry)) {
852  result = tidal_power(deriv_zone_index, false, entry);
853  } else {
854  for(unsigned i = 0; i < __orbit_entries.size(); ++i)
855  if(entry == __orbit_entries[i])
856  return __orbit_power[i];
857 
858  assert(false);
859  }
860  if(__above_lock_fractions.size() == 0)
861  return result;
862  } else {
863  double above_frac = (__above_lock_fractions
865  [zone(deriv_zone_index).locked_zone_index()]);
866  result = (above_frac
867  *
868  tidal_power(deriv_zone_index, true, entry)
869  +
870  (1.0 - above_frac)
871  *
872  tidal_power(deriv_zone_index, false, entry));
873  }
874 
875  assert(above_lock_fraction_deriv.size()
876  ==
878 
879  unsigned correction_index = 0;
880  for(unsigned zone_index = 0; zone_index < number_zones(); ++zone_index)
881  if(zone(zone_index).locked()) {
882  result += (
883  above_lock_fraction_deriv
884  [zone(zone_index).locked_zone_index()]
885  *
886  __orbit_power_correction[correction_index]
887  );
888  ++correction_index;
889  }
890  return result;
891  }
892 
894  Dissipation::QuantityEntry entry,
895  unsigned deriv_zone_index,
896  const Eigen::VectorXd &above_lock_fraction_deriv
897  ) const
898  {
899  if(zone_specific(entry) && deriv_zone_index != 0) {
900  const DissipatingZone &deriv_zone = zone(deriv_zone_index);
901  double above_frac = (__above_lock_fractions
903  [deriv_zone.locked_zone_index()]);
904  Eigen::Vector3d result = zone_to_zone_transform(
905  deriv_zone,
906  zone(0),
907  (
908  (above_frac - 1.0)
909  *
910  __tidal_torques_below[deriv_zone_index][entry]
911  -
912  above_frac
913  *
914  __tidal_torques_above[deriv_zone_index][entry]
915  )
916  );
917 
918  assert(above_lock_fraction_deriv.size()
919  ==
921 
922  unsigned correction_index = 0;
923  for(
924  unsigned zone_index = 0;
925  zone_index < number_zones();
926  ++zone_index
927  )
928  if(zone(zone_index).locked()) {
929  result += (
930  above_lock_fraction_deriv[
931  zone(zone_index).locked_zone_index()
932  ]
933  *
934  __orbit_torque_correction[correction_index]
935  );
936  ++correction_index;
937  }
938  if(
939  entry == Dissipation::INCLINATION
940  ||
941  entry == Dissipation::PERIAPSIS
942  )
943  result += zone_to_zone_transform(
944  deriv_zone,
945  zone(0),
946  (
947  above_frac
948  *
949  (
951  [deriv_zone_index]
953  )
954  +
955  (1.0 - above_frac)
956  *
957  (
959  [deriv_zone_index]
960  [Dissipation::NO_DERIV]
961  )
962  ),
963  entry,
964  true
965  );
966  return result;
967  } else {
968  return __orbit_torque[entry];
969  }
970  }
971 
973  const DissipatingZone &reference_zone,
974  Dissipation::QuantityEntry entry,
975  unsigned deriv_zone_index,
976  const Eigen::VectorXd &above_lock_fraction_deriv
977  ) const
978  {
979  Eigen::Vector3d result;
980  result = zone_to_zone_transform(
981  zone(0),
982  reference_zone,
983  tidal_orbit_torque(entry,
984  deriv_zone_index,
985  above_lock_fraction_deriv)
986  );
987 
988  if(
989  (
990  entry == Dissipation::INCLINATION
991  ||
992  entry == Dissipation::PERIAPSIS
993  )
994  &&
995  (
996  deriv_zone_index == 0
997  ||
998  &zone(deriv_zone_index) == &reference_zone
999  )
1000  ) {
1001  result += zone_to_zone_transform(zone(0),
1002  reference_zone,
1004  entry,
1005  deriv_zone_index==0);
1006  }
1007  return result;
1008  }
1009 
1011  {
1012  for(unsigned zone_ind = 0; zone_ind < number_zones(); ++zone_ind)
1013  zone(zone_ind).add_to_evolution();
1014  }
1015 
1017  {
1018  for(unsigned zone_ind = 0; zone_ind < number_zones(); ++zone_ind)
1019  zone(zone_ind).rewind_evolution(nsteps);
1020  }
1021 
1023  {
1024  for(unsigned zone_ind = 0; zone_ind < number_zones(); ++zone_ind)
1025  zone(zone_ind).reset_evolution();
1026  }
1027 
1029  BinarySystem &system,
1030  bool primary
1031  )
1032  {
1034  for(unsigned zone_ind = 0; zone_ind < number_zones(); ++zone_ind)
1035  (*result) |= zone(zone_ind).stopping_conditions(system,
1036  primary,
1037  zone_ind);
1038  return result;
1039  }
1040 
1041  void DissipatingBody::change_expansion_order(unsigned new_expansion_order,
1042  BinarySystem &system,
1043  bool primary)
1044  {
1045  for(unsigned zone_ind = 0; zone_ind < number_zones(); ++zone_ind)
1046  if(zone(zone_ind).dissipative())
1047  zone(zone_ind).change_expansion_order(new_expansion_order,
1048  system,
1049  primary,
1050  zone_ind);
1051  }
1052 
1053 } //End Envolve namespace.
double __power_norm
The coefficient used to normalize tidal power.
double tidal_power(bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const
The dimensionless tidal power or one of its derivatives.
unsigned locked_zone_index() const
The index of this zone in the list of locked zones (valid only if locked).
std::valarray< std::valarray< Eigen::Vector3d > > __angular_momentum_transfer
The rate of angular momentum transfer between two neighboring zones due to zone boundaries moving...
virtual void rewind_evolution(unsigned nsteps)
Discards the last steps from the evolution.
virtual void reset_evolution()
Discards all evolution.
RADIUS
The derivative w.r.t. the radius of the body in .
virtual void configure(bool initialize, double age, double orbital_frequency, double eccentricity, double orbital_angmom, double spin, double inclination, double periapsis, bool spin_is_frequency)
Defines the current orbit, triggering re-calculation of all quantities.
virtual void add_to_evolution()
Appends the state defined by last configure(), to the evolution.
double angular_momentum() const
The angular momentum of the given zone in .
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.
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...
SPIN_FREQUENCY
The derivative w.r.t. the spin frequency of a dissipating zone.
DissipatingBody()
Some initializations for new objects.
double tidal_power(unsigned zone_index, bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const
Tidal power dissipated in the given zone.
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...
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.
virtual void add_to_evolution()
Appends the state defined by last configure(), to the evolution.
virtual const DissipatingZone & zone(unsigned zone_index) const =0
A modifiable reference to one of the body&#39;s zones.
SEMIMAJOR
The derivative w.r.t. the semimajor axis in AU.
MOMENT_OF_INERTIA
Eigen::Vector3d zone_to_zone_transform(const ZoneOrientation &from_zone, const ZoneOrientation &to_zone, const Eigen::Vector3d &vector, Dissipation::QuantityEntry deriv, bool with_respect_to_from)
Transforms a vector betwen the coordinates systems of two zones.
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 tidal_torque_y(bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const
The dimensionless torque along y.
NUM_DERIVATIVES
The total number of derivatives supported.
virtual double moment_of_inertia(int deriv_order=0) const =0
Moment of inertia of the zone or its age derivative at the age of last configure() call...
double mass() const
The mass of the body (constant with age).
void calculate_orbit_rate_corrections()
Calculates the correction the orbital energy gain and torque due to switching locked zones from fully...
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.
virtual CombinedStoppingCondition * stopping_conditions(BinarySystem &system, bool primary)
Conditions detecting the next possible discontinuities in the evolution due to this body...
double __dorbital_frequency_da
The derivative of the orbital frequency w.r.t. semimajor axis.
double radius(int deriv_order=0) const
The current radius or its derivative with age of the body.
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()).
const double solar_radius
Radius of the sun [m].
A layer of a system body for which the tidal bulge is not exactly in phase with the tidal potential...
virtual void reset_evolution()
Discards all evolution.
double __orbital_frequency
The mean angular velocity with which the orbit is traversed.
double normalize_torques(double companion_mass, double semimajor, double orbital_frequency)
Scales the dimensionless torques as appropriate and corrects the relevant derivatives returning the n...
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...
const double Gyr
Gyr [s].
virtual double outer_mass(int deriv_order=0) const =0
Mass coordinate of the zone&#39;s outer ouboundary or its derivative.
void correct_orbit_torque(std::valarray< Eigen::VectorXd > &above_lock_fractions)
Corrects __orbit_torque for zone locks.
ECCENTRICITY
The derivative w.r.t. the eccentricity.
END_DIMENSIONLESS_DERIV
const double day
day [s]
const double solar_mass
Mass of the sun [kg].
virtual void change_expansion_order(unsigned new_expansion_order, BinarySystem &system, bool primary, unsigned zone_index)
Changes the order of the eccentricity expansion performed.
AGE
The derivative w.r.t. age, excluding the dependence through the body&#39;s radius and the moments of iner...
std::vector< Eigen::Vector3d > __orbit_torque_correction
Corrections to __orbit_torque_below (undifferentiated) if single zones switch to above.
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.
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.
virtual CombinedStoppingCondition * stopping_conditions(BinarySystem &system, bool primary, unsigned zone_index)
Conditions detecting the next possible discontinuities in the evolution due to this zone...
virtual bool locked(int orbital_frequency_multiplier, int spin_frequency_multiplier) const
Should return true iff the given term is presently locked.
NO_DERIV
The quantity itself, undifferentiated.
Declares the DissipatingBody class.
ORBITAL_FREQUENCY
The derivative w.r.t. the orbital frequency.
double tidal_torque_z(bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const
The dimensionless tidal torque along z.
std::vector< Dissipation::QuantityEntry > __orbit_entries
The quantities w.r.t. which derivatives of the orbit energy gain and torque are pre-calculated.
SPIN_ANGMOM
The derivative w.r.t. the spin angular momentum in .
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()).
virtual void rewind_evolution(unsigned nsteps)
Discards the last steps from the evolution.
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()).
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...
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&#39;s tidal dis...
std::valarray< std::valarray< Eigen::Vector3d > > __tidal_torques_above
Tidal torques and their derivatives on each zone for spin frequency approaching potential lock from a...
A class combining the the outputs of multiple stopping conditions.
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.
std::vector< Eigen::Vector3d > __orbit_torque
Torque on the orbit (and its derivatives and error) due to tides on this body.
virtual bool dissipative() const =0
Should return true iff the zone has some non-zero dissipation.
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.
double tidal_torque_x(bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const
The dimensionless tidal torque along x.
Describes a system of two bodies orbiting each other.
Definition: BinarySystem.h:57
std::valarray< double > __orbit_power
Total tidal power (and derivatives and error) gained by the orbit from this body. ...
virtual void change_expansion_order(unsigned new_expansion_order, BinarySystem &system, bool primary)
Change the eccentricity expansion order for all dissipative zones.
std::valarray< double > __orbit_power_correction
Corrections to __orbit_power_below (undifferentiated) if single zones switch to above.
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.
virtual double outer_radius(int deriv_order=0) const =0
Outer radius of the zone or its derivative (per last.
double spin_frequency() const
The surface spin freuqency of the body.
double spin_frequency() const
The spin frequency of the given zone.
const double G
Gravitational constant in SI.
virtual unsigned number_zones() const =0
The number of zones the body consists of.