Planetary Orbital Evolution due to Tides
Orbital evolution of two objects experiencing tides
BinarySystem.cpp
Go to the documentation of this file.
1 
8 #define BUILDING_LIBRARY
9 #include "BinarySystem.h"
10 
11 namespace Evolve {
12 
14  {
15  __locked_zones.clear();
16  for(short body_ind = 0; body_ind < 2; ++body_ind) {
17  DissipatingBody &body = (body_ind == 0 ? __body1 : __body2);
18  for(
19  unsigned zone_ind = 0;
20  zone_ind < body.number_zones();
21  ++zone_ind
22  )
23  if(body.zone(zone_ind).locked()) {
24  body.zone(zone_ind).locked_zone_index() =
25  __locked_zones.size();
26  __locked_zones.push_back(zone_ind);
27  }
28  }
29  }
30 
31  BinarySystem::lock_scenario_type BinarySystem::find_synchronized_zones(
32  double precision
33  )
34  {
35  lock_scenario_type result;
36  double worb = orbital_frequency();
37  for(unsigned short body_ind = 0; body_ind < 2; ++body_ind) {
38  DissipatingBody &body = (body_ind == 0 ? __body1 : __body2);
39  for(
40  unsigned zone_ind = 0;
41  zone_ind < body.number_zones();
42  ++zone_ind
43  ) {
44  DissipatingZone &zone = body.zone(zone_ind);
45  if(!zone.can_lock())
46  continue;
47  double wspin = zone.spin_frequency();
48  bool other_lock=false;
49  unsigned scenario_zone_i = (zone_ind
50  +
51  body_ind * __body1.number_zones());
52  do {
53  const SpinOrbitLockInfo &lock = zone.lock_monitored(
54  other_lock
55  );
56  if(
57  std::abs(
58  (
59  lock.orbital_frequency_multiplier() * worb
60  -
61  lock.spin_frequency_multiplier() * wspin
62  )
63  /
64  worb
65  ) < precision
66  ) {
67  if(other_lock)
68  zone.swap_monitored_locks();
69  result.push_back(
70  std::make_tuple(
71  scenario_zone_i,
72  lock.lock_direction(),
73  zone.angular_momentum()
74  )
75  );
76  } else
77  assert(!zone.locked()
78  ||
79  std::get<0>(result.back()) == scenario_zone_i);
80  other_lock = !other_lock;
81  } while(other_lock);
82  }
83  }
84  return result;
85  }
86 
88  const lock_scenario_type &lock_scenario
89  )
90  {
91  assert(__evolution_mode == Core::BINARY);
92 #ifdef VERBOSE_DEBUG
93  std::cerr << "Setting lock scenario: " << std::endl;
94  describe_lock_scenario(std::cerr,
95  lock_scenario,
96  std::vector<bool>(lock_scenario.size(), false),
97  true);
98 #endif
99 
100  std::vector<double>
101  angmom(number_zones()),
102  inclinations(number_zones()),
103  periapses(number_zones() - 1);
104 
105  std::vector<double>::iterator
106  angmom_dest = angmom.begin(),
107  inclination_dest = inclinations.begin(),
108  periapsis_dest = periapses.begin();
109  lock_scenario_type::const_iterator
110  scenario_zone_iter = lock_scenario.begin();
111 
112  unsigned scenario_zone;
113  short scenario_lock_dir;
114  double scenario_angmom, zone_angmom;
115  std::tie(scenario_zone,
116  scenario_lock_dir,
117  scenario_angmom) = *scenario_zone_iter;
118 
119  for(short body_i = 0; body_i < 2; ++body_i) {
120  DissipatingBody &body = (body_i == 0 ? __body1 : __body2);
121  for(unsigned zone_i = 0; zone_i < body.number_zones(); ++zone_i) {
122  DissipatingZone &zone = body.zone(zone_i);
123  const SpinOrbitLockInfo &current_lock = zone.lock_monitored();
124  if(
125  zone_i + body_i * __body1.number_zones()
126  ==
127  scenario_zone
128  ) {
129  zone_angmom = scenario_angmom;
130  if(
131  current_lock.lock_direction()
132  !=
133  scenario_lock_dir
134  ) {
135  if(zone.locked()) {
136  body.unlock_zone_spin(zone_i, scenario_lock_dir);
137 #ifdef VERBOSE_DEBUG
138  std::cerr << "Unlocking "
139  << (body_i == 0 ? "primary" : "secondary")
140  << " zone " << zone_i
141  << "(overall zone " << scenario_zone
142  << ")"
143  << std::endl;
144 #endif
145  } else {
146  body.lock_zone_spin(
147  zone_i,
148  current_lock.orbital_frequency_multiplier(),
149  current_lock.spin_frequency_multiplier()
150  );
151 #ifdef VERBOSE_DEBUG
152  std::cerr << "Locking "
153  << (body_i == 0 ? "primary" : "secondary")
154  << " zone " << zone_i
155  << "(overall zone " << scenario_zone
156  << ")"
157  << std::endl;
158 #endif
159 
160  if(scenario_lock_dir != 0) {
161 #ifdef VERBOSE_DEBUG
162  std::cerr << "Unlocking "
163  << (body_i == 0 ? "primary" : "secondary")
164  << " zone " << zone_i
165  << "(overall zone " << scenario_zone
166  << ")"
167  << std::endl;
168 #endif
169  body.unlock_zone_spin(zone_i,
170  scenario_lock_dir);
171  }
172  }
173  }
174 #ifdef VERBOSE_DEBUG
175  else {
176  std::cerr << (body_i == 0 ? "Primary" : "Secondary")
177  << " zone " << zone_i
178  << "(overall zone " << scenario_zone
179  << ") already matches scenario."
180  << std::endl;
181  }
182 #endif
183  ++scenario_zone_iter;
184  if(scenario_zone_iter != lock_scenario.end())
185  std::tie(scenario_zone,
186  scenario_lock_dir,
187  scenario_angmom) = *scenario_zone_iter;
188 
189  } else {
190  zone_angmom = zone.angular_momentum();
191  }
192 
193  if(!zone.locked())
194  *angmom_dest++ = zone_angmom;
195  *inclination_dest++ = zone.inclination();
196  if(body_i || zone_i)
197  *periapsis_dest++ = zone.periapsis();
198  }
199  }
200  configure(false,
201  __age,
202  __semimajor,
204  angmom.data(),
205  inclinations.data(),
206  periapses.data(),
207  Core::BINARY);
208  }
209 
211  const std::vector<short> &unlock_directions,
212  const std::vector<double> &original_angmom
213  )
214  {
215  std::vector<double>
216  config_angmom(number_zones()),
217  config_inclinations(number_zones()),
218  config_periapses(number_zones() - 1);
219  unsigned locked_i = 0;
220  for(unsigned config_i = 0; config_i < number_zones(); ++config_i) {
221  DissipatingBody &body = (config_i < __body1.number_zones()
222  ? __body1
223  : __body2);
224  unsigned body_zone_i = (
225  config_i < __body1.number_zones()
226  ? config_i
227  : config_i - __body1.number_zones()
228  );
229  DissipatingZone &zone = body.zone(body_zone_i);
230 
231  config_inclinations[config_i] = zone.inclination();
232  if(config_i)
233  config_periapses[config_i - 1] = zone.periapsis();
234 
235  if(zone.locked()) {
236  config_angmom[config_i] = (original_angmom.size()
237  ? original_angmom[locked_i]
238  : zone.angular_momentum());
239  body.unlock_zone_spin(body_zone_i,
240  unlock_directions[locked_i]);
241  ++locked_i;
242  } else
243  config_angmom[config_i] = zone.angular_momentum();
244  }
245  configure(false,
246  __age,
247  __semimajor,
249  config_angmom.data(),
250  config_inclinations.data(),
251  config_periapses.data(),
252  Core::BINARY);
253  }
254 
256  {
257  bool primary_zone = test_zone_ind < __body1.number_zones();
258  unsigned body_zone_ind = (primary_zone
259  ? test_zone_ind
260  : test_zone_ind - __body1.number_zones());
261  const SpinOrbitLockInfo &lock = (
262  primary_zone
263  ? __body1
264  : __body2
265  ).zone(
266  body_zone_ind
267  ).lock_monitored();
268  std::valarray<double> orbit;
269  fill_orbit(orbit);
270  std::valarray<double> derivatives(orbit.size());
272  &(orbit[0]),
273  Core::BINARY,
274  &(derivatives[0]));
275 
276  std::valarray<double> check_cond_deriv;
278  lock,
279  primary_zone,
280  body_zone_ind,
281  *this
282  )(
283  Core::BINARY,
284  orbit,
285  derivatives,
286  check_cond_deriv
287  );
288  assert(check_cond_deriv.size() == 1);
289 #ifndef NDEBUG
290  std::cerr << "Lock dir " << lock.lock_direction()
291  << "for zone " << test_zone_ind
292  << " (condition deriv: " << check_cond_deriv[0]
293  << (lock.lock_direction() * check_cond_deriv[0] > 0
294  ? ") is good"
295  : ") is bad")
296  << std::endl;
297 #endif
298  return lock.lock_direction() * check_cond_deriv[0] < 0;
299  }
300 
301 #ifndef NDEBUG
303  std::ostream &os,
304  const lock_scenario_type &lock_scenario,
305  const std::vector<bool> passed,
306  bool add_header
307  )
308  {
309  lock_scenario_type::const_iterator
310  scenario_zone_i = lock_scenario.begin();
311  std::vector<bool>::const_iterator passed_i = passed.begin();
312  unsigned test_zone_ind;
313  short lock_direction;
314  std::tie(test_zone_ind, lock_direction, std::ignore) = *scenario_zone_i;
315 
316  if(add_header)
317  os << "||" << std::setw(9 * __body1.number_zones() - 1) << "Body1"
318  << "||" << std::setw(9 * __body2.number_zones() - 1) << "Body2"
319  << "||" << std::endl
320  << std::string(
321  9 * (__body1.number_zones() + __body2.number_zones()) + 4,
322  '_'
323  )
324  << std::endl
325  << "||";
326 
327  for(short body_i = 0; body_i < 2; ++body_i) {
328  DissipatingBody &body = (body_i == 0 ? __body1 : __body2);
329  for(unsigned zone_i = 0; zone_i < body.number_zones(); ++zone_i) {
330  if(zone_i != 0)
331  os << "|";
332  DissipatingZone &zone = body.zone(zone_i);
333  os << " ";
334  if(
335  zone_i + body_i * __body1.number_zones()
336  ==
337  test_zone_ind
338  ) {
339  if(lock_direction == 0)
340  os << "0";
341  else if(lock_direction > 0)
342  os << "+";
343  os << lock_direction
344  << " (" << (*passed_i ? "V" : "X") << ")";
345  ++scenario_zone_i;
346  if(scenario_zone_i != lock_scenario.end())
347  std::tie(test_zone_ind,
348  lock_direction,
349  std::ignore) = *scenario_zone_i;
350  } else {
351  if(zone.can_lock()) os << " NSZ ";
352  else os << " NLZ ";
353  }
354  os << " ";
355  }
356  os << "||";
357  }
358  os << std::endl;
359  }
360 #endif
361 
363  const lock_scenario_type &lock_scenario
364 #ifndef NDEBUG
365  , bool first_scenario
366 #endif
367  )
368  {
369  set_lock_scenario(lock_scenario);
370  double *above_lock_fraction_p = __above_lock_fractions[
372  ].data();
373 #ifndef DEBUG
374  std::vector<bool> passed(lock_scenario.size());
375  std::vector<bool>::iterator passed_i = passed.begin();
376  bool result = true;
377 #endif
378  for(
379  lock_scenario_type::const_iterator
380  scenario_zone_iter = lock_scenario.begin();
381  scenario_zone_iter != lock_scenario.end();
382  ++scenario_zone_iter
383  ) {
384  unsigned test_zone_ind;
385  short lock_direction;
386  std::tie(test_zone_ind,
387  lock_direction,
388  std::ignore) = *scenario_zone_iter;
389  if(lock_direction == 0) {
390 #ifdef VERBOSE_DEBUG
391  std::cerr << "Lock for zone " << test_zone_ind << " will ";
392 #endif
393  if(
394  !(
395  *above_lock_fraction_p > 0
396  &&
397  *above_lock_fraction_p < 1
398  )
399  ) {
400 #ifndef NDEBUG
401  *passed_i++ = false;
402  result = false;
403 #else
404  return false;
405 #endif
406 #ifdef VERBOSE_DEBUG
407  std::cerr << "not ";
408 #endif
409  }
410 #ifndef NDEBUG
411  else *passed_i++ = true;
412 #endif
413 #ifdef VERBOSE_DEBUG
414  std::cerr << "hold (above lock frac: "
415  << *above_lock_fraction_p
416  << ")" << std::endl;
417 #endif
418  ++above_lock_fraction_p;
419  } else if(!test_synchronized_unlocked_zone(test_zone_ind)) {
420 #ifndef NDEBUG
421  *passed_i++ = false;
422  result = false;
423 #else
424  return false;
425 #endif
426  }
427 #ifndef NDEBUG
428  else *passed_i++ = true;
429 #endif
430  }
431 #ifndef NDEBUG
432  describe_lock_scenario(std::cerr,
433  lock_scenario,
434  passed,
435  first_scenario || true);
436  return result;
437 #else
438  return true;
439 #endif
440  }
441 
443  lock_scenario_type::const_iterator next_synchronized_zone,
444  unsigned num_synchronized_zones,
445  lock_scenario_type &lock_scenario
446 #ifndef NDEBUG
447  , bool first_scenario
448 #endif
449 
450  )
451  {
452  if(lock_scenario.size() == num_synchronized_zones) {
453 #ifndef NDEBUG
454  if(test_lock_scenario(lock_scenario, first_scenario)) {
455  __selected_lock_scenario = lock_scenario;
456  std::cerr << "Selected lock scenario: " << std::endl;
458  std::cerr,
459  lock_scenario,
460  std::vector<bool>(lock_scenario.size(), true),
461  true
462  );
464  std::cerr,
466  std::vector<bool>(lock_scenario.size(), true),
467  true
468  );
469 
470  return true;
471  } else
472  return false;
473 #else
474  return test_lock_scenario(lock_scenario);
475 #endif
476  }
477 
478  assert(lock_scenario.size() < num_synchronized_zones);
479 
480  lock_scenario.push_back(*next_synchronized_zone);
481  short &scenario_dir = std::get<1>(lock_scenario.back());
482  short orig_lock_direction = std::get<1>(*next_synchronized_zone);
483 #ifndef NDEBUG
484  bool result = false;
485 #endif
486  ++next_synchronized_zone;
487  for(scenario_dir = -1; scenario_dir <= 1; ++scenario_dir) {
488  if(
490  next_synchronized_zone,
491  num_synchronized_zones,
492  lock_scenario
493 #ifndef NDEBUG
494  , first_scenario
495 #endif
496  )
497  ) {
498 #ifndef NDEBUG
499  assert(!result);
500  result = true;
501 #else
502  return true;
503 #endif
504  }
505 #ifndef NDEBUG
506  first_scenario = false;
507 #endif
508  }
509  lock_scenario.pop_back();
510 #ifndef NDEBUG
511  return result;
512 #else
513  return false;
514 #endif
515  }
516 
517 
519  double *evolution_rates
520  ) const
521  {
522  DissipatingZone &locked_zone = __body1.zone(0);
523  locked_zone.set_evolution_rates(
524  locked_zone.moment_of_inertia(1) * locked_zone.spin_frequency(),
525  0.0,
526  0.0
527  );
528  for(
529  unsigned zone_index = 1;
530  zone_index < __body1.number_zones();
531  ++zone_index
532  ) {
533  evolution_rates[zone_index - 1] = __body1.nontidal_torque(
534  zone_index
535  )[2];
536 
537  __body1.zone(zone_index).set_evolution_rates(
538  evolution_rates[zone_index - 1],
539  0.0,
540  0.0
541  );
542 
543  assert(__body1.nontidal_torque(zone_index)[0] == 0);
544  assert(__body1.nontidal_torque(zone_index)[1] == 0);
545 
546  }
547  return 0;
548  }
549 
550 #ifdef ENABLE_DERIVATIVES
551  void BinarySystem::locked_surface_jacobian(double *param_derivs,
552  double *age_derivs) const
553  {
554  unsigned num_param = __body1.number_zones() - 1;
555  double dR_dt = __body1.radius(1),
556  dIabove_dt = __body1.zone(0).moment_of_inertia(1),
557  dI_dt = __body1.zone(1).moment_of_inertia(1);
558  for(unsigned row = 0; row < num_param; ++row) {
559  double dIbelow_dt = 0;
560  age_derivs[row] =
562  +
563  __body1.nontidal_torque(row + 1, Dissipation::RADIUS)[2] * dR_dt
564  +
565  __body1.nontidal_torque(row + 1,
567  0)[2] * dI_dt
568  +
569  __body1.nontidal_torque(row + 1,
571  -1)[2] * dIabove_dt;
572  if(row < num_param - 1) {
573  dIbelow_dt = __body1.zone(row + 2).moment_of_inertia(1);
574  age_derivs[row] += __body1.nontidal_torque(
575  row + 1,
577  1)[2] * dIbelow_dt;
578  }
579  dIabove_dt = dI_dt;
580  dI_dt = dIbelow_dt;
581 
582  assert(__body1.nontidal_torque(row + 1, Dissipation::AGE)[0] == 0);
583  assert(__body1.nontidal_torque(row + 1, Dissipation::AGE)[1] == 0);
584  assert(__body1.nontidal_torque(row + 1, Dissipation::RADIUS)[0]
585  ==
586  0);
587  assert(__body1.nontidal_torque(row + 1, Dissipation::RADIUS)[1]
588  ==
589  0);
590  assert(__body1.nontidal_torque(row + 1,
592  -1)[0]
593  ==
594  0);
595  assert(__body1.nontidal_torque(row + 1,
597  -1)[1]
598  ==
599  0);
600  assert(__body1.nontidal_torque(row + 1,
602  0)[0]
603  ==
604  0);
605  assert(__body1.nontidal_torque(row + 1,
607  0)[1]
608  ==
609  0);
610  assert(__body1.nontidal_torque(row + 1,
612  1)[0]
613  ==
614  0);
615  assert(__body1.nontidal_torque(row + 1,
617  1)[1]
618  ==
619  0);
620 
621  for(unsigned col = 0; col < num_param; ++col) {
622  double &dest = *(param_derivs + row * num_param + col);
623  if(std::abs(static_cast<int>(row) - static_cast<int>(col)) < 1)
624  dest = 0;
625  else {
626  dest = __body1.nontidal_torque(row + 1,
628  col - row)[2];
629 
630  assert(__body1.nontidal_torque(row + 1,
632  col - row)[0]
633  ==
634  0);
635  assert(__body1.nontidal_torque(row + 1,
637  col - row)[1]
638  ==
639  0);
640 
641  }
642  }
643  }
644  }
645 #endif
646 
648  double *evolution_rates
649  ) const
650  {
651  unsigned nzones = __body1.number_zones();
652  double ref_angmom = __body1.zone(0).angular_momentum(),
653  *inclination_evol = evolution_rates,
654  *periapsis_evol = evolution_rates + (nzones - 1),
655  *angmom_evol = periapsis_evol+(nzones - 1);
656 
657  Eigen::Vector3d reference_torque;
658  for(unsigned zone_index = 0; zone_index < nzones; ++zone_index) {
659  Eigen::Vector3d torque = __body1.nontidal_torque(zone_index);
660  angmom_evol[zone_index] = torque[2];
661 
662  DissipatingZone &zone = __body1.zone(zone_index);
663  if(zone_index) {
664  zone.set_reference_zone_angmom(ref_angmom);
665  inclination_evol[zone_index - 1] = zone.inclination_evolution(
667  zone,
668  reference_torque),
669  torque
670  );
671  periapsis_evol[zone_index - 1] = zone.periapsis_evolution(
673  zone,
674  reference_torque),
675  torque
676  );
677  } else reference_torque = torque;
678 
679  zone.set_evolution_rates(
680  angmom_evol[zone_index],
681  (zone_index ? inclination_evol[zone_index - 1] : 0.0),
682  (zone_index ? periapsis_evol[zone_index -1] : 0.0)
683  );
684  }
685 #ifdef VERBOSE_DEBUG
686  std::cerr << "rates: ";
687  for(unsigned i = 0; i < 3 * nzones - 2; ++i) {
688  if(i) std::cerr << ", ";
689  std::cerr << evolution_rates[i];
690  }
691  std::cerr << std::endl;
692 #endif
693  return 0;
694  }
695 
696 #ifdef ENABLE_DERIVATIVES
697  void BinarySystem::fill_single_body_jacobian(
698  double *inclination_param_derivs,
699  double *periapsis_param_derivs,
700  double *angmom_param_derivs,
701  double *inclination_age_derivs,
702  double *periapsis_age_derivs,
703  double *angmom_age_derivs
704  ) const
705  {
706  unsigned nzones=__body1.number_zones(), nparams = 3 * nzones - 2;
707  Eigen::Vector3d ref_torque = __body1.nontidal_torque(0),
708  dref_torque_dincl = __body1.nontidal_torque(
709  0,
711  1
712  ),
713  dref_torque_dperi=__body1.nontidal_torque(
714  0,
716  1
717  ),
718  dref_torque_dangmom0=__body1.nontidal_torque(
719  0,
721  0
722  ),
723  dref_torque_dangmom1=__body1.nontidal_torque(
724  0,
726  1
727  ),
728  dref_torque_dage = __body1.nontidal_torque(
729  0,
731  ),
732  zero3d(0, 0, 0);
733  for(
734  unsigned deriv_zone_ind = 0;
735  deriv_zone_ind < nzones;
736  ++deriv_zone_ind
737  ) {
738  angmom_param_derivs[deriv_zone_ind - 1] =
739  (deriv_zone_ind == 1 ? dref_torque_dincl[2] : 0);
740  angmom_param_derivs[nzones+deriv_zone_ind - 2] =
741  (deriv_zone_ind == 1 ? dref_torque_dperi[2] : 0);
742  double &angmom_angmom_deriv =
743  angmom_param_derivs[2 * nzones+deriv_zone_ind - 2];
744  if(deriv_zone_ind==0) angmom_angmom_deriv = dref_torque_dangmom0[2];
745  else if(deriv_zone_ind == 1)
746  angmom_angmom_deriv = dref_torque_dangmom1[2];
747  else angmom_angmom_deriv = 0;
748  }
749  angmom_age_derivs[0] = dref_torque_dage[2];
750  DissipatingZone &ref_zone = __body1.zone(0);
751  for(unsigned zone_ind = 1; zone_ind < nzones; ++zone_ind) {
752  DissipatingZone &zone = __body1.zone(zone_ind);
753  Eigen::Vector3d
754  zone_ref_torque = zone_to_zone_transform(ref_zone,
755  zone,
756  ref_torque),
757  zone_torque = __body1.nontidal_torque(zone_ind),
758  ref_torque_age_deriv = zone_to_zone_transform(ref_zone,
759  zone,
760  dref_torque_dage),
761  zone_torque_age_deriv=__body1.nontidal_torque(zone_ind,
763  0);
764  inclination_age_derivs[zone_ind - 1] =
765  zone.inclination_evolution(zone_ref_torque,
766  zone_torque,
768  ref_torque_age_deriv,
769  zone_torque_age_deriv);
770  periapsis_age_derivs[zone_ind - 1] =
771  zone.periapsis_evolution(zone_ref_torque,
772  zone_torque,
774  ref_torque_age_deriv,
775  zone_torque_age_deriv);
776  angmom_age_derivs[zone_ind] = zone_torque_age_deriv[2];
777  for(
778  unsigned quantity_ind = 0;
779  quantity_ind < nparams;
780  ++quantity_ind
781  ) {
782  unsigned offset = (zone_ind - 1) * nparams + quantity_ind;
783  double &inclination_dest = inclination_param_derivs[offset],
784  &periapsis_dest = periapsis_param_derivs[offset],
785  &angmom_dest = angmom_param_derivs[offset+nparams];
786  unsigned quantity_zone = quantity_ind;
787  Dissipation::QuantityEntry with_respect_to;
788  Eigen::Vector3d ref_torque_deriv(0, 0, 0);
789  if(quantity_ind > 2 * nzones - 2) {
790  quantity_zone -= 2 * nzones - 2;
791  with_respect_to = Dissipation::SPIN_ANGMOM;
792  if(quantity_zone == 0)
793  ref_torque_deriv = dref_torque_dangmom0;
794  else if(quantity_zone == 1)
795  ref_torque_deriv = dref_torque_dangmom1;
796  if(quantity_zone <= 1)
797  ref_torque_deriv = zone_to_zone_transform(
798  ref_zone,
799  zone,
800  ref_torque_deriv
801  );
802  } else {
803  if(quantity_ind >= nzones - 1) {
804  quantity_zone -= nzones - 2;
805  with_respect_to = Dissipation::PERIAPSIS;
806  if(quantity_zone == 1)
807  ref_torque_deriv = dref_torque_dperi;
808  } else {
809  quantity_zone += 1;
810  with_respect_to = Dissipation::INCLINATION;
811  if(quantity_zone == 1)
812  ref_torque_deriv = dref_torque_dincl;
813  }
814  if(quantity_zone == 1)
815  ref_torque_deriv = zone_to_zone_transform(
816  ref_zone,
817  zone,
818  ref_torque_deriv
819  );
820  ref_torque_deriv += zone_to_zone_transform(ref_zone,
821  zone,
822  ref_torque,
823  with_respect_to);
824  }
825  Eigen::Vector3d zone_torque_deriv;
826  if(quantity_zone == zone_ind) {
827  zone_torque_deriv = __body1.nontidal_torque(zone_ind,
828  with_respect_to,
829  0);
830  inclination_dest=zone.inclination_evolution(
831  zone_ref_torque,
832  zone_torque,
833  with_respect_to,
834  ref_torque_deriv,
835  zone_torque_deriv
836  );
837  periapsis_dest=zone.periapsis_evolution(zone_ref_torque,
838  zone_torque,
839  with_respect_to,
840  ref_torque_deriv,
841  zone_torque_deriv);
842  } else {
843  if(std::abs(static_cast<int>(quantity_zone-zone_ind)) == 1) {
844  zone_torque_deriv = __body1.nontidal_torque(
845  zone_ind,
846  with_respect_to,
847  quantity_zone - zone_ind
848  );
849  } else zone_torque_deriv = zero3d;
850  inclination_dest = zone.inclination_evolution(
851  ref_torque_deriv,
852  zone_torque_deriv
853  );
854  periapsis_dest = zone.periapsis_evolution(ref_torque_deriv,
855  zone_torque_deriv);
856  }
857  angmom_dest = zone_torque_deriv[2];
858  }
859  }
860  }
861 
862  void BinarySystem::single_body_jacobian(double *param_derivs,
863  double *age_derivs) const
864  {
865  unsigned nzones = __body1.number_zones(),
866  nparams=3*nzones-2;
867  fill_single_body_jacobian(param_derivs,
868  param_derivs + (nzones - 1) * nparams,
869  param_derivs + 2 * (nzones - 1) * nparams,
870  age_derivs,
871  age_derivs + (nzones - 1),
872  age_derivs + 2 * (nzones - 1));
873  }
874 #endif
875 
877  double orbit_power,
878  double orbit_power_deriv
879  ) const
880  {
881  if(std::isnan(orbit_power_deriv))
882  return (orbit_power == 0
883  ? 0
884  : -__semimajor * orbit_power / __orbital_energy);
885 
886  else if(orbit_power == 0 && orbit_power_deriv == 0)
887  return 0;
888 
889  return (
890  -(2.0 * orbit_power
891  +
892  __semimajor * orbit_power_deriv
893  )
894  /
896  );
897  }
898 
899  double BinarySystem::eccentricity_evolution(double orbit_power,
900  double orbit_angmom_gain,
901  double orbit_power_deriv,
902  double orbit_angmom_gain_deriv,
903  bool semimajor_deriv) const
904  {
905  if(__eccentricity <= 1e-8) return 0;
906  double e2 = std::pow(__eccentricity, 2),
907  factor = -(1.0 - e2) / (2.0 * __eccentricity);
908 
909  if(std::isnan(orbit_power_deriv))
910  return factor * (orbit_power / __orbital_energy
911  +
912  2.0 * orbit_angmom_gain / __orbital_angmom);
913  else if(semimajor_deriv)
914  return (
915  factor
916  *
917  (
918  (
919  orbit_power_deriv
920  +
921  orbit_power / __semimajor
922  )
923  /
925  +
926  (
927  2.0 * orbit_angmom_gain_deriv
928  -
929  orbit_angmom_gain / __semimajor
930  )
931  /
933  )
934  );
935  else return (
936  factor
937  *
938  (
939  orbit_power_deriv / __orbital_energy
940  +
941  2.0 * orbit_angmom_gain_deriv / __orbital_angmom
942  )
943  -
944  2.0 * orbit_angmom_gain / __orbital_angmom
945  -
946  (1.0 + e2) / e2 * (
947  orbit_power / __orbital_energy
948  +
949  2.0 * orbit_angmom_gain / __orbital_angmom
950  )
951  );
952  }
953 
955  Dissipation::QuantityEntry entry,
956  bool body1_deriv,
957  Eigen::MatrixXd &matrix,
958  Eigen::VectorXd &rhs
959  ) const
960  {
961  unsigned deriv_zone_index = (body1_deriv
962  ? 0
964  DissipatingBody &deriv_body = (body1_deriv ? __body1 : __body2);
965  DissipatingZone &deriv_zone = deriv_body.zone(0);
966 
967  if(
969  ||
970  entry == Dissipation::SEMIMAJOR
971  ) {
972  double coef = (entry == Dissipation::SEMIMAJOR
973  ? 1
974  : Core::orbital_angular_velocity(__body1.mass(),
975  __body2.mass(),
976  __semimajor,
977  true));
978  bool done = false;
979  unsigned locked_ind = 0;
980  for(DissipatingBody *body = &__body1; !done; body = &__body2) {
981  for(
982  unsigned zone_ind = 0;
983  zone_ind < body->number_zones();
984  ++zone_ind
985  ) {
986  rhs.array() += (1.5
987  *
988  (
990  +
992  )
993  /
995  /
997  *
998  coef);
999  if(body->zone(zone_ind).locked()) {
1000  matrix.col(locked_ind).array() += (
1001  1.5
1002  *
1003  (
1004  body->tidal_power(zone_ind, true)
1005  -
1006  body->tidal_power(zone_ind, false)
1007  )
1008  /
1010  /
1011  __semimajor * coef
1012  );
1013 
1014  ++locked_ind;
1015  }
1016  }
1017  done = (body == &__body2);
1018  }
1019  } else if(deriv_zone.locked()) {
1020  if(
1022  ||
1023  entry == Dissipation::SPIN_ANGMOM
1024  ) {
1025  double coef = (
1027  ? deriv_zone.moment_of_inertia()
1028  : 1
1029  ) / std::pow(deriv_zone.angular_momentum(), 2);
1030  matrix(deriv_zone_index, deriv_zone_index) -=(
1031  deriv_body.tidal_torque(0, true)[2]
1032  -
1033  deriv_body.tidal_torque(0, false)[2]
1034  ) * coef;
1035  rhs(deriv_zone_index) -= (
1036  deriv_body.tidal_torque(0, false)[2]
1037  *
1038  coef
1039  );
1040  } else if(entry == Dissipation::MOMENT_OF_INERTIA) {
1041  rhs(deriv_zone_index) -= (
1042  deriv_zone.moment_of_inertia(1)
1043  /
1044  std::pow(deriv_zone.moment_of_inertia(), 2)
1045  );
1046  }
1047  }
1048  }
1049 
1051  Eigen::VectorXd &fractions,
1052  Dissipation::QuantityEntry entry,
1053  bool body1_deriv
1054  )
1055  {
1056  unsigned num_locked_zones = __locked_zones.size();
1057 
1058  assert(num_locked_zones == (__body1.number_locked_zones()
1059  +
1061 #ifndef NDEBUG
1062  std::cerr << "Setting "
1063  << num_locked_zones
1064  << " above lock fractions "
1065  << entry
1066  << std::endl;
1067 #endif
1068 
1069  if(num_locked_zones == 0) {
1070  fractions.resize(0);
1071  return;
1072  }
1073  std::valarray<double> nontidal_torque(num_locked_zones),
1074  tidal_torque_z_above(num_locked_zones),
1075  tidal_torque_z_below(num_locked_zones),
1076  tidal_power_difference(num_locked_zones);
1077  unsigned locked_zone_ind = 0;
1078  for(
1079  std::list<unsigned>::const_iterator zi = __locked_zones.begin();
1080  zi != __locked_zones.end();
1081  ++zi
1082  ) {
1083  DissipatingBody &body = (
1084  locked_zone_ind < __body1.number_locked_zones()
1085  ? __body1
1086  : __body2
1087  );
1088  nontidal_torque[locked_zone_ind] =
1089  body.nontidal_torque(*zi, entry)[2];
1090  tidal_torque_z_above[locked_zone_ind] =
1091  body.tidal_torque(*zi, true, entry)[2];
1092  tidal_torque_z_below[locked_zone_ind] =
1093  body.tidal_torque(*zi, false, entry)[2];
1094  if(!zone_specific(entry) || *zi == 0) {
1095  tidal_power_difference[locked_zone_ind] = (
1096  body.tidal_power(*zi, true, entry)
1097  -
1098  body.tidal_power(*zi, false, entry)
1099  );
1100  }
1101  ++locked_zone_ind;
1102  }
1103  Eigen::MatrixXd matrix(num_locked_zones, num_locked_zones);
1104  Eigen::VectorXd rhs(num_locked_zones);
1105  rhs.setConstant(1.5
1106  *
1107  (
1108  __body1.tidal_orbit_power(entry)
1109  +
1110  __body2.tidal_orbit_power(entry)
1111  )
1112  /
1114  unsigned i = 0;
1115  for(
1116  std::list<unsigned>::const_iterator
1117  zi=__locked_zones.begin();
1118  zi!=__locked_zones.end();
1119  ++zi
1120  ) {
1121  if(!zone_specific(entry) || *zi == 0)
1122  matrix.col(i).setConstant(1.5
1123  *
1124  tidal_power_difference[i]
1125  /
1127  else matrix.col(i).setZero();
1129  ? __body1
1130  : __body2).zone(*zi);
1131  matrix(i, i) += (
1132  (tidal_torque_z_above[i] - tidal_torque_z_below[i])
1133  /
1134  zone.angular_momentum()
1135  );
1136  if(entry == Dissipation::NO_DERIV)
1137  rhs(i) += (zone.moment_of_inertia(1)
1138  /
1139  zone.moment_of_inertia());
1140  else if(entry == Dissipation::AGE)
1141  rhs(i) += (zone.moment_of_inertia(2)
1142  /
1143  zone.moment_of_inertia());
1144  rhs(i) -= ((nontidal_torque[i] + tidal_torque_z_below[i])
1145  /
1146  zone.angular_momentum());
1147  ++i;
1148  }
1149  above_lock_problem_deriv_correction(entry, body1_deriv, matrix, rhs);
1150  if(entry == Dissipation::NO_DERIV) {
1151  __above_lock_fractions_decomp.compute(matrix);
1152  fractions = __above_lock_fractions_decomp.solve(rhs);
1153  } else {
1154  fractions = __above_lock_fractions_decomp.solve(
1156  );
1157  }
1158  }
1159 
1161  Dissipation::QuantityEntry entry,
1162  DissipatingBody &body,
1163  unsigned zone_index)
1164  {
1165  assert(entry == Dissipation::INCLINATION
1166  ||
1167  entry == Dissipation::PERIAPSIS
1168  ||
1170  ||
1171  entry == Dissipation::SPIN_ANGMOM);
1172  assert(zone_index > 0 || &body == &__body2);
1173  assert(body.number_zones() > zone_index);
1174  assert(number_locked_zones()
1175  ==
1177 
1178  unsigned num_locked_zones = number_locked_zones();
1179  if(num_locked_zones == 0) return Eigen::VectorXd();
1180  DissipatingZone &deriv_zone = body.zone(zone_index);
1181  unsigned locked_zone_index=(&body == &__body1
1182  ? 0
1184  Eigen::VectorXd rhs;
1185  for(unsigned i = 0; i < zone_index; ++i)
1186  if(body.zone(i).locked()) ++locked_zone_index;
1187  if(deriv_zone.locked()) {
1188  double above_frac =
1189  __above_lock_fractions[Dissipation::NO_DERIV][locked_zone_index];
1190  rhs.setConstant(
1191  num_locked_zones,
1192  -1.5 * (above_frac*body.tidal_power(zone_index, true, entry)
1193  +
1194  (1.0 - above_frac) * body.tidal_power(zone_index,
1195  false,
1196  entry))
1197  -
1198  body.nontidal_torque(zone_index, entry)[2]
1199  );
1200  rhs(locked_zone_index) -=
1201  above_frac*body.tidal_torque(zone_index, true, entry)[2]
1202  +
1203  (1.0 - above_frac) * body.tidal_torque(zone_index,
1204  false,
1205  entry)[2];
1206  if(entry == Dissipation::MOMENT_OF_INERTIA)
1207  rhs(locked_zone_index) -= (
1208  deriv_zone.moment_of_inertia(1)
1209  /
1210  std::pow(deriv_zone.moment_of_inertia(), 2)
1211  );
1212  else if(entry == Dissipation::SPIN_ANGMOM)
1213  rhs(locked_zone_index) += (
1214  above_frac*body.tidal_torque(zone_index, true)[2]
1215  +
1216  (1.0 - above_frac)
1217  *
1218  body.tidal_torque(zone_index, false)[2]
1219  +
1220  body.nontidal_torque(zone_index)[2]
1221  ) / std::pow(deriv_zone.angular_momentum(), 2);
1222  } else {
1223  rhs.setConstant(num_locked_zones,
1224  -1.5 * body.tidal_power(zone_index, false, entry));
1225  }
1226  if(zone_index > 0 && body.zone(zone_index - 1).locked())
1227  rhs(body.zone(zone_index - 1).locked_zone_index()) -=
1228  body.nontidal_torque(zone_index - 1, entry, 1)[2];
1229  if(
1230  zone_index + 1 < body.number_zones()
1231  &&
1232  body.zone(zone_index + 1).locked()
1233  )
1234  rhs(body.zone(zone_index + 1).locked_zone_index()) -=
1235  body.nontidal_torque(zone_index + 1, entry, -1)[2];
1236  return __above_lock_fractions_decomp.solve(rhs);
1237  }
1238 
1240  {
1241  unsigned num_zones = __body1.number_zones() + __body2.number_zones();
1242  DissipatingBody *body = &__body1;
1243  __above_lock_fractions_inclination_deriv.resize(num_zones);
1244  __above_lock_fractions_periapsis_deriv.resize(num_zones);
1245  __above_lock_fractions_inertia_deriv.resize(num_zones);
1246  __above_lock_fractions_angmom_deriv.resize(num_zones);
1255  unsigned body_zone_ind = 0;
1256  for(unsigned zone_ind = 1; zone_ind < num_zones; ++zone_ind) {
1257  ++body_zone_ind;
1258  if(body_zone_ind == body->number_zones()) {
1259  body_zone_ind = 0;
1260  body = &__body2;
1261  }
1264  *body,
1265  body_zone_ind);
1268  *body,
1269  body_zone_ind);
1272  *body,
1273  body_zone_ind);
1276  *body,
1277  body_zone_ind);
1278  }
1279  }
1280 
1282  {
1283  std::valarray<Eigen::VectorXd>
1284  body2_above_lock_fractions(Dissipation::NUM_DERIVATIVES);
1285  for(
1286  int deriv = Dissipation::NO_DERIV;
1288  ++deriv
1289  ) {
1291  __above_lock_fractions[deriv],
1292  static_cast<Dissipation::QuantityEntry>(deriv)
1293  );
1294  if(!zone_specific(static_cast < Dissipation::QuantityEntry>(deriv)))
1295  body2_above_lock_fractions[deriv] = __above_lock_fractions[deriv];
1297  body2_above_lock_fractions[deriv],
1298  static_cast<Dissipation::QuantityEntry>(deriv),
1299  false
1300  );
1301  }
1303  body2_above_lock_fractions[Dissipation::RADIUS];
1305  __body2.set_above_lock_fractions(body2_above_lock_fractions);
1307  }
1308 
1310  {
1311  __orbit_torque = (
1313  +
1315  __body1.zone(0),
1317  );
1318 
1320  +
1322 
1324  *
1325  std::sin(__body1.zone(0).inclination())
1326  +
1327  __orbit_torque[2]
1328  *
1329  std::cos(__body1.zone(0).inclination()));
1330  }
1331 
1333  double *differential_equations
1334  ) const
1335  {
1336  DissipatingZone &reference_zone = __body1.zone(0);
1337  unsigned num_body1_zones = __body1.number_zones(),
1338  num_total_zones = num_body1_zones + __body2.number_zones();
1339  double *inclination_rates = differential_equations + 2,
1340  *periapsis_rates = inclination_rates + num_total_zones,
1341  *angmom_rates = periapsis_rates + num_total_zones - 1;
1342  unsigned angmom_skipped = 0;
1343  double reference_periapsis_rate = Core::NaN;
1344  for(unsigned zone_ind = 0; zone_ind < num_total_zones; ++zone_ind) {
1345  DissipatingBody &body = (zone_ind < num_body1_zones
1346  ? __body1
1347  : __body2);
1348  unsigned body_zone_ind = zone_ind;
1349  if(zone_ind >= num_body1_zones)
1350  body_zone_ind -= num_body1_zones;
1351  DissipatingZone &zone = body.zone(body_zone_ind);
1352  Eigen::Vector3d zone_orbit_torque,
1353  total_zone_torque;
1354  total_zone_torque = body.nontidal_torque(body_zone_ind);
1355 
1356  zone_orbit_torque = (zone_ind
1357  ? zone_to_zone_transform(reference_zone,
1358  zone,
1360  : __orbit_torque);
1361 
1362  Dissipation::QuantityEntry entry = Dissipation::NO_DERIV;
1363  if(zone.locked()) {
1364  double above_frac = (__above_lock_fractions
1366  [angmom_skipped]);
1367  total_zone_torque += (
1368  above_frac
1369  *
1370  body.tidal_torque(body_zone_ind, true, entry)
1371  +
1372  (1.0 - above_frac)
1373  *
1374  body.tidal_torque(body_zone_ind, false, entry)
1375  );
1376  ++angmom_skipped;
1377  } else total_zone_torque += body.tidal_torque(body_zone_ind,
1378  false,
1379  entry);
1380  inclination_rates[zone_ind] = zone.inclination_evolution(
1381  zone_orbit_torque,
1382  total_zone_torque,
1383  entry
1384  );
1385  assert(!std::isnan(inclination_rates[zone_ind]));
1386  if(zone_ind) {
1387  periapsis_rates[zone_ind - 1] = zone.periapsis_evolution(
1388  zone_orbit_torque,
1389  total_zone_torque,
1390  entry
1391  ) - reference_periapsis_rate;
1392  assert(!std::isnan(periapsis_rates[zone_ind - 1]));
1393  } else {
1394  reference_periapsis_rate = zone.periapsis_evolution(
1395  zone_orbit_torque,
1396  total_zone_torque,
1397  entry
1398  );
1399  assert(!std::isnan(reference_periapsis_rate));
1400 
1401  }
1402  if(!zone.locked()) {
1403  angmom_rates[zone_ind - angmom_skipped] = total_zone_torque[2];
1404  assert(!std::isnan(angmom_rates[zone_ind - angmom_skipped]));
1405  }
1406 
1407  zone.set_evolution_rates(
1408  total_zone_torque[2],
1409  inclination_rates[zone_ind],
1410  (zone_ind ? periapsis_rates[zone_ind - 1] : 0.0)
1411  );
1412 
1413 #ifdef VERBOSE_DEBUG
1414  std::cerr << "Zone " << zone_ind
1415  << " torque: " << total_zone_torque
1416  << " inclination rate";
1417  std::cerr << ": " << inclination_rates[zone_ind];
1418  if(zone_ind)
1419  std::cerr << " periapsis rate: "
1420  << periapsis_rates[zone_ind - 1];
1421 
1422  std::cerr << std::endl;
1423 #endif
1424  assert(!std::isnan(total_zone_torque.sum()));
1425 
1426  }
1427 
1428  differential_equations[0] = semimajor_evolution(__orbit_power);
1429  differential_equations[1] = eccentricity_evolution(__orbit_power,
1431  __semimajor_rate = differential_equations[0];
1432  __eccentricity_rate = differential_equations[1];
1433 
1434  if(!angmom_skipped)
1435  differential_equations[0] *= 6.5 * std::pow(__semimajor, 5.5);
1436 
1437 #ifdef VERBOSE_DEBUG
1438  std::cerr << "rates: ";
1439  for(
1440  unsigned i = 0;
1441  i < 3 * (__body1.number_zones() + __body2.number_zones()) - 1;
1442  ++i
1443  ) {
1444  if(i) std::cerr << ", ";
1445  std::cerr << differential_equations[i];
1446  }
1447  std::cerr << std::endl;
1448 #endif
1449 
1450  return 0;
1451  }
1452 
1453  template<typename VALUE_TYPE>
1455  const DissipatingBody &body,
1456  VALUE_TYPE (DissipatingBody::*func)(Dissipation::QuantityEntry,
1457  unsigned,
1458  const Eigen::VectorXd &) const,
1459  std::valarray<VALUE_TYPE> &orbit_rate_deriv,
1460  unsigned offset
1461  ) const
1462  {
1463  unsigned num_zones = __body1.number_zones() + __body2.number_zones(),
1464  num_param = orbit_rate_deriv.size() - 1;
1465  orbit_rate_deriv[0] += (body.*func)(Dissipation::SEMIMAJOR,
1466  0,
1467  Eigen::VectorXd());
1468  orbit_rate_deriv[1] += (body.*func)(Dissipation::ECCENTRICITY,
1469  0,
1470  Eigen::VectorXd());
1471  orbit_rate_deriv[num_param] += (
1472  (body.*func)(Dissipation::AGE, 0, Eigen::VectorXd())
1473  +
1474  (body.*func)(Dissipation::RADIUS, 0, Eigen::VectorXd())
1475  *
1476  body.radius(1)
1477  );
1478  unsigned locked_zone_count = (offset ? __body1.number_locked_zones() : 0);
1479  for(unsigned zone_ind = 0; zone_ind < body.number_zones(); ++zone_ind) {
1480  orbit_rate_deriv[zone_ind+2+offset] =
1481  (body.*func)(
1483  zone_ind,
1484  __above_lock_fractions_inclination_deriv[zone_ind + offset]
1485  );
1486  if(zone_ind+offset > 0)
1487  orbit_rate_deriv[zone_ind + 1 + num_zones + offset] =
1488  (body.*func)(
1490  zone_ind,
1492  );
1493  const DissipatingZone &zone = body.zone(zone_ind);
1494  if(zone.locked()) ++locked_zone_count;
1495  else orbit_rate_deriv[
1496  zone_ind + 1 + 2 * num_zones - locked_zone_count
1497  ] = (body.*func)(
1499  zone_ind,
1500  __above_lock_fractions_angmom_deriv[zone_ind + offset]
1501  );
1502  orbit_rate_deriv[num_param] += (
1503  (body.*func)(
1505  zone_ind,
1506  __above_lock_fractions_inertia_deriv[zone_ind + offset]
1507  )
1508  *
1509  zone.moment_of_inertia(1)
1510  );
1511  }
1512  }
1513 
1515  std::valarray<double> &orbit_power_deriv) const
1516  {
1517  orbit_power_deriv[0] = 0;
1518  orbit_power_deriv[1] = 0;
1519  orbit_power_deriv[orbit_power_deriv.size() - 1] = 0;
1521  orbit_power_deriv,
1522  0);
1525  orbit_power_deriv,
1526  __body1.number_zones());
1527  }
1528 
1530  std::valarray<double> &orbit_angmom_gain_deriv
1531  ) const
1532  {
1533  std::valarray<Eigen::Vector3d>
1534  body1_orbit_torque_deriv(orbit_angmom_gain_deriv.size()),
1535  body2_orbit_torque_deriv(orbit_angmom_gain_deriv.size());
1536  body1_orbit_torque_deriv[0] =
1537  body1_orbit_torque_deriv[1] =
1538  body2_orbit_torque_deriv[0] =
1539  body2_orbit_torque_deriv[1] = Eigen::Vector3d(0, 0, 0);
1542  body1_orbit_torque_deriv,
1543  0);
1545  body2_orbit_torque_deriv,
1546  __body1.number_zones());
1547  double body1_sin_inc = std::sin(__body1.zone(0).inclination()),
1548  body1_cos_inc = std::sin(__body1.zone(0).inclination()),
1549  body2_sin_inc = std::sin(__body2.zone(0).inclination()),
1550  body2_cos_inc = std::sin(__body2.zone(0).inclination());
1551  for(unsigned i = 0; i < orbit_angmom_gain_deriv.size(); ++i)
1552  orbit_angmom_gain_deriv[i] = (
1553  body1_sin_inc * body1_orbit_torque_deriv[i][0]
1554  +
1555  body1_cos_inc * body1_orbit_torque_deriv[i][2]
1556  +
1557  body2_sin_inc * body2_orbit_torque_deriv[i][0]
1558  +
1559  body2_cos_inc * body2_orbit_torque_deriv[i][2]
1560  );
1561  Eigen::Vector3d body1_orbit_torque = __body1.tidal_orbit_torque(),
1562  body2_orbit_torque = __body2.tidal_orbit_torque();
1563  orbit_angmom_gain_deriv[2] +=(body1_cos_inc * body1_orbit_torque[0]
1564  -
1565  body1_sin_inc * body1_orbit_torque[2]
1566  +
1567  body2_cos_inc * body2_orbit_torque[0]
1568  -
1569  body2_sin_inc * body2_orbit_torque[2]);
1570  }
1571 
1572 #ifdef ENABLE_DERIVATIVES
1573  void BinarySystem::semimajor_jacobian(
1574  const std::valarray<double> &orbit_power_deriv,
1575  bool a6p5,
1576  double *param_derivs,
1577  double &age_deriv
1578  ) const
1579  {
1580  param_derivs[0] = semimajor_evolution(__orbit_power,
1581  orbit_power_deriv[0]);
1582  if(a6p5) param_derivs[0] +=
1584  unsigned i = 1;
1585  for(; i < orbit_power_deriv.size() - 1; ++i)
1586  param_derivs[i] = semimajor_evolution(orbit_power_deriv[i]);
1587  age_deriv = semimajor_evolution(orbit_power_deriv[i]);
1588  }
1589 
1590  void BinarySystem::eccentricity_jacobian(
1591  const std::valarray<double> &orbit_power_deriv,
1592  const std::valarray<double> &orbit_angmom_gain_deriv,
1593  bool a6p5,
1594  double *param_derivs,
1595  double &age_deriv
1596  ) const
1597  {
1598  param_derivs[0] = eccentricity_evolution(__orbit_power,
1600  orbit_power_deriv[0],
1601  orbit_angmom_gain_deriv[0],
1602  true);
1603  if(a6p5) param_derivs[0] /= 6.5 * std::pow(__semimajor, 5.5);
1604  param_derivs[1] = eccentricity_evolution(__orbit_power,
1606  orbit_power_deriv[1],
1607  orbit_angmom_gain_deriv[1],
1608  false);
1609  unsigned i = 2;
1610  for(; i < orbit_power_deriv.size() - 1; ++i)
1611  param_derivs[i] = eccentricity_evolution(orbit_power_deriv[i],
1612  orbit_angmom_gain_deriv[i]);
1613  age_deriv=eccentricity_evolution(orbit_power_deriv[i],
1614  orbit_angmom_gain_deriv[i]);
1615  }
1616 
1617  void BinarySystem::angle_evolution_age_deriv(DissipatingBody &body,
1618  unsigned zone_ind,
1619  double sin_inc,
1620  double cos_inc,
1621  unsigned locked_zone_ind,
1622  double &inclination,
1623  double &periapsis) const
1624  {
1625  double dR1_dt = __body1.radius(1),
1626  dR2_dt = __body2.radius(1);
1627  DissipatingZone &zone = body.zone(zone_ind);
1628  Eigen::Vector3d orbit_torque_age_deriv = (
1630  +
1632  +
1634  +
1636  );
1637  DissipatingBody *other_body = &__body1;
1638  unsigned body_zone_ind = 0,
1639  num_zones = number_zones();
1640  for(
1641  unsigned other_zone_ind = 0;
1642  other_zone_ind<num_zones;
1643  ++other_zone_ind
1644  ) {
1645  orbit_torque_age_deriv += (
1646  other_body->tidal_orbit_torque(
1647  zone,
1649  body_zone_ind,
1650  __above_lock_fractions_inertia_deriv[other_zone_ind]
1651  )
1652  *
1653  other_body->zone(body_zone_ind).moment_of_inertia(1)
1654  );
1655  ++body_zone_ind;
1656  if(body_zone_ind == other_body->number_zones()) {
1657  body_zone_ind = 0;
1658  other_body = &__body2;
1659  }
1660  }
1661  periapsis = (orbit_torque_age_deriv[1]
1662  *
1663  cos_inc
1664  /
1665  (__orbital_angmom * sin_inc));
1666  inclination = (orbit_torque_age_deriv[0] * cos_inc
1667  -
1668  orbit_torque_age_deriv[2] * sin_inc) / __orbital_angmom;
1669  Eigen::Vector3d to_add = -body.nontidal_torque(zone_ind,
1671  if(zone.locked()) {
1672  double above_frac =
1674  double above_frac_deriv =
1675  __above_lock_fractions[Dissipation::AGE][locked_zone_ind];
1676  to_add -= (
1677  body.tidal_torque(zone_ind, true, Dissipation::AGE)
1678  *
1679  above_frac
1680  +
1681  body.tidal_torque(zone_ind, false, Dissipation::AGE)
1682  *
1683  (1.0 - above_frac)
1684  +
1685  above_frac_deriv*(body.tidal_torque(zone_ind, true)
1686  -
1687  body.tidal_torque(zone_ind, false))
1688  );
1689  } else to_add += body.tidal_torque(zone_ind, false, Dissipation::AGE);
1690  to_add /= zone.angular_momentum();
1691  inclination += to_add[0];
1692  periapsis -= to_add[1] / sin_inc;
1693  }
1694 
1695  void BinarySystem::angle_evolution_orbit_deriv(Dissipation::QuantityEntry entry,
1696  double angmom_deriv,
1697  DissipatingBody &body,
1698  unsigned zone_ind,
1699  double sin_inc,
1700  double cos_inc,
1701  unsigned locked_zone_ind,
1702  double &inclination,
1703  double &periapsis) const
1704  {
1705  assert(entry == Dissipation::SEMIMAJOR
1706  ||
1707  entry == Dissipation::ECCENTRICITY);
1708 
1709  DissipatingZone &zone = body.zone(zone_ind);
1710  Eigen::Vector3d orbit_torque = (__body1.tidal_orbit_torque(zone)
1711  +
1712  __body2.tidal_orbit_torque(zone)),
1713  orbit_torque_deriv = (
1714  __body1.tidal_orbit_torque(zone, entry)
1715  +
1716  __body2.tidal_orbit_torque(zone, entry)
1717  );
1718  inclination = (orbit_torque_deriv[0] * cos_inc
1719  -
1720  orbit_torque_deriv[2] * sin_inc
1721  -
1722  (orbit_torque[0] * cos_inc - orbit_torque[2] * sin_inc)
1723  *
1724  angmom_deriv
1725  /
1727  periapsis = (
1728  (
1729  orbit_torque[1] * angmom_deriv / __orbital_angmom
1730  -
1731  orbit_torque_deriv[1]
1732  )
1733  *
1734  cos_inc
1735  /
1736  (__orbital_angmom * sin_inc)
1737  );
1738  Eigen::Vector3d to_add =- body.nontidal_torque(zone_ind, entry);
1739  if(zone.locked()) {
1740  double
1741  above_frac =
1743  above_frac_deriv = __above_lock_fractions[entry][locked_zone_ind];
1744  to_add -= (
1745  body.tidal_torque(zone_ind, true, entry) * above_frac
1746  +
1747  body.tidal_torque(zone_ind, false, entry)
1748  *
1749  (1.0 - above_frac)
1750  +
1751  above_frac_deriv
1752  *
1753  (
1754  body.tidal_torque(zone_ind, true)
1755  -
1756  body.tidal_torque(zone_ind, false)
1757  )
1758  );
1759  } else to_add -= body.tidal_torque(zone_ind, false, entry);
1760  to_add /= zone.angular_momentum();
1761  inclination += to_add[0];
1762  periapsis -= to_add[1] / sin_inc;
1763  }
1764 
1765  void BinarySystem::fill_orbit_torque_deriv(
1766  Dissipation::QuantityEntry entry,
1767  DissipatingBody &body,
1768  unsigned zone_ind,
1769  std::valarray<Eigen::Vector3d> &orbit_torque_deriv
1770  ) const
1771  {
1772  assert(entry == Dissipation::INCLINATION
1773  ||
1774  entry == Dissipation::PERIAPSIS
1775  ||
1777  ||
1778  entry == Dissipation::SPIN_ANGMOM);
1779  assert(orbit_torque_deriv.size()
1780  ==
1782 
1783  unsigned body_deriv_zone_ind = 0,
1784  num_zones = orbit_torque_deriv.size();
1785  DissipatingBody *deriv_body = &__body1;
1786  const std::valarray<Eigen::VectorXd> *above_frac_deriv;
1787  if(entry == Dissipation::INCLINATION)
1788  above_frac_deriv = &__above_lock_fractions_inclination_deriv;
1789  else if(entry == Dissipation::PERIAPSIS)
1790  above_frac_deriv = &__above_lock_fractions_periapsis_deriv;
1791  else if(entry == Dissipation::MOMENT_OF_INERTIA)
1792  above_frac_deriv = &__above_lock_fractions_inertia_deriv;
1793  else
1794  above_frac_deriv = &__above_lock_fractions_angmom_deriv;
1795 
1796  DissipatingZone &zone = body.zone(zone_ind);
1797  for(
1798  unsigned deriv_zone_ind = 0;
1799  deriv_zone_ind < num_zones;
1800  ++deriv_zone_ind
1801  ) {
1802  orbit_torque_deriv[deriv_zone_ind] =
1803  deriv_body->tidal_orbit_torque(
1804  zone,
1805  entry,
1806  body_deriv_zone_ind,
1807  (*above_frac_deriv)[deriv_zone_ind]
1808  );
1809  if(
1810  (
1811  entry == Dissipation::INCLINATION
1812  ||
1813  entry == Dissipation::PERIAPSIS
1814  )
1815  &&
1816  deriv_body == &body
1817  &&
1818  zone_ind == body_deriv_zone_ind
1819  ) {
1820  DissipatingBody &other_body = (&body == &__body1
1821  ? __body2
1822  : __body1);
1823  orbit_torque_deriv[deriv_zone_ind] +=
1824  zone_to_zone_transform(other_body.zone(0),
1825  zone,
1826  other_body.tidal_orbit_torque(),
1827  entry,
1828  false);
1829  }
1830  ++body_deriv_zone_ind;
1831  if(body_deriv_zone_ind == deriv_body->number_zones()) {
1832  body_deriv_zone_ind = 0;
1833  deriv_body = &__body2;
1834  }
1835  }
1836  }
1837 
1838  void BinarySystem::fill_zone_torque_deriv(
1839  Dissipation::QuantityEntry entry,
1840  DissipatingBody &body,
1841  unsigned zone_ind,
1842  std::valarray<Eigen::Vector3d> &zone_torque_deriv
1843  ) const
1844  {
1845 #ifndef NDEBUG
1846  if(body.zone(zone_ind).locked())
1847  assert(zone_torque_deriv.size() == 4);
1848  else
1849  assert(zone_torque_deriv.size() == 3);
1850 #endif
1851 
1852  zone_torque_deriv[0] = (zone_ind == 0
1853  ? Eigen::Vector3d::Zero()
1854  : body.nontidal_torque(zone_ind, entry, -1));
1855  Eigen::Vector3d nontidal_torque = body.nontidal_torque(zone_ind,
1856  entry,
1857  0);
1858  zone_torque_deriv[1] = (nontidal_torque
1859  +
1860  body.tidal_torque(zone_ind, false, entry));
1861  zone_torque_deriv[2] = (zone_ind < body.number_zones() - 1
1862  ? body.nontidal_torque(zone_ind, entry, 1)
1863  : Eigen::Vector3d::Zero());
1864  if(body.zone(zone_ind).locked())
1865  zone_torque_deriv[3] = (nontidal_torque
1866  +
1867  body.tidal_torque(zone_ind, true, entry));
1868  }
1869 
1870  void BinarySystem::inclination_evolution_zone_derivs(
1871  Dissipation::QuantityEntry entry,
1872  DissipatingBody &body,
1873  unsigned zone_ind,
1874  double zone_x_torque_above,
1875  double zone_x_torque_below,
1876  const std::valarray<Eigen::Vector3d> &zone_torque_deriv,
1877  const Eigen::Vector3d &orbit_torque,
1878  const std::valarray<Eigen::Vector3d> &orbit_torque_deriv,
1879  const std::valarray<Eigen::VectorXd> &above_frac_deriv,
1880  double sin_inc,
1881  double cos_inc,
1882  unsigned locked_zone_ind,
1883  double *result
1884  ) const
1885  {
1886  assert(entry == Dissipation::INCLINATION
1887  ||
1888  entry == Dissipation::PERIAPSIS
1889  ||
1890  entry == Dissipation::SPIN_ANGMOM);
1891  assert(orbit_torque_deriv.size()
1892  ==
1894 
1895  DissipatingZone &zone = body.zone(zone_ind);
1896  unsigned num_zones = orbit_torque_deriv.size();
1897  unsigned global_zone_ind = zone_ind;
1898  if(&body == &__body2) global_zone_ind += __body1.number_zones();
1899 
1900  double above_frac =
1902 
1903  for(
1904  unsigned deriv_zone_ind = (entry == Dissipation::PERIAPSIS ? 1 : 0);
1905  deriv_zone_ind < num_zones;
1906  ++deriv_zone_ind
1907  ) {
1908  double &dest = (entry == Dissipation::PERIAPSIS
1909  ? result[deriv_zone_ind - 1]
1910  : result[deriv_zone_ind]);
1911  dest = (
1912  orbit_torque_deriv[deriv_zone_ind][0] * cos_inc
1913  -
1914  orbit_torque_deriv[deriv_zone_ind][2] * sin_inc
1915  ) / __orbital_angmom;
1916  if(zone.locked())
1917  dest -= (above_frac_deriv[deriv_zone_ind][locked_zone_ind]
1918  *
1919  (zone_x_torque_above - zone_x_torque_below)
1920  /
1921  zone.angular_momentum());
1922  if(std::abs(static_cast<int>(deriv_zone_ind - global_zone_ind)) <= 1) {
1923  double torque_deriv=
1924  zone_torque_deriv[deriv_zone_ind + 1 - global_zone_ind][0];
1925  if(zone.locked() && deriv_zone_ind == global_zone_ind) {
1926  dest-=(
1927  above_frac * zone_torque_deriv[3][0]
1928  +
1929  (1.0 - above_frac) * torque_deriv
1930  ) / zone.angular_momentum();
1931  } else dest -= torque_deriv / zone.angular_momentum();
1932  }
1933  }
1934  if(entry == Dissipation::INCLINATION)
1935  result[global_zone_ind] -= (
1936  orbit_torque[0] * sin_inc
1937  +
1938  orbit_torque[2] * cos_inc
1939  ) / __orbital_angmom;
1940  else if(entry == Dissipation::SPIN_ANGMOM)
1941  result[global_zone_ind] += (
1942  zone.locked()
1943  ? (above_frac * zone_x_torque_above
1944  +
1945  (1.0 - above_frac) * zone_x_torque_below)
1946  : zone_x_torque_below
1947  ) / std::pow(zone.angular_momentum(), 2);
1948  }
1949 
1950  void BinarySystem::periapsis_evolution_zone_derivs(
1951  Dissipation::QuantityEntry entry,
1952  DissipatingBody &body,
1953  unsigned zone_ind,
1954  double zone_y_torque_above,
1955  double zone_y_torque_below,
1956  const std::valarray<Eigen::Vector3d> &zone_torque_deriv,
1957  double orbit_y_torque,
1958  const std::valarray<Eigen::Vector3d> &orbit_torque_deriv,
1959  const std::valarray<Eigen::VectorXd> &above_frac_deriv,
1960  double sin_inc,
1961  double cos_inc,
1962  unsigned locked_zone_ind,
1963  double *result
1964  ) const
1965  {
1966  assert(entry == Dissipation::INCLINATION
1967  ||
1968  entry == Dissipation::PERIAPSIS
1969  ||
1971  ||
1972  entry == Dissipation::SPIN_ANGMOM);
1973  assert(orbit_torque_deriv.size()
1974  ==
1976 
1977  DissipatingZone &zone=body.zone(zone_ind);
1978  unsigned num_zones = orbit_torque_deriv.size();
1979  unsigned global_zone_ind = zone_ind;
1980  if(&body == &__body2) global_zone_ind += __body1.number_zones();
1981 
1982  double above_frac =
1984 
1985  for(
1986  unsigned deriv_zone_ind = 0;
1987  deriv_zone_ind < num_zones;
1988  ++deriv_zone_ind
1989  ) {
1990  if(entry == Dissipation::PERIAPSIS && deriv_zone_ind == 0) continue;
1991  double &dest = (entry == Dissipation::PERIAPSIS
1992  ? result[deriv_zone_ind - 1]
1993  : result[deriv_zone_ind]);
1994  dest = (-orbit_torque_deriv[deriv_zone_ind][1]
1995  *
1996  cos_inc
1997  /
1998  (__orbital_angmom * sin_inc));
1999  if(zone.locked())
2000  dest += (above_frac_deriv[deriv_zone_ind][locked_zone_ind]
2001  *
2002  (zone_y_torque_above - zone_y_torque_below)
2003  /
2004  (sin_inc * zone.angular_momentum()));
2005  if(
2006  std::abs(static_cast<int>(deriv_zone_ind - global_zone_ind))
2007  <=
2008  1
2009  ) {
2010  double torque_deriv =
2011  zone_torque_deriv[deriv_zone_ind + 1 - global_zone_ind][1];
2012  if(zone.locked() && deriv_zone_ind == global_zone_ind) {
2013  dest+=(
2014  above_frac * zone_torque_deriv[3][1]
2015  +
2016  (1.0 - above_frac) * torque_deriv
2017  ) / (sin_inc * zone.angular_momentum());
2018  } else dest += (torque_deriv
2019  /
2020  (sin_inc * zone.angular_momentum()));
2021  }
2022  }
2023  double zone_torque;
2024  if(zone.locked()) zone_torque = (above_frac*zone_y_torque_above
2025  +
2026  (1.0 - above_frac)
2027  *
2028  zone_y_torque_below);
2029  else zone_torque = zone_y_torque_below;
2030  if(entry == Dissipation::INCLINATION)
2031  result[global_zone_ind] += (
2032  orbit_y_torque / (__orbital_angmom * std::pow(sin_inc, 2))
2033  -
2034  zone_torque*cos_inc
2035  /
2036  (zone.angular_momentum() * std::pow(sin_inc, 2))
2037  );
2038  else if(entry == Dissipation::SPIN_ANGMOM)
2039  result[global_zone_ind] -=
2040  zone_torque / (std::pow(zone.angular_momentum(), 2) * sin_inc);
2041  }
2042 
2043  void BinarySystem::spin_angmom_evolution_zone_derivs(
2044  Dissipation::QuantityEntry entry,
2045  DissipatingBody &body,
2046  unsigned zone_ind,
2047  double zone_z_torque_above,
2048  double zone_z_torque_below,
2049  const std::valarray<Eigen::Vector3d> &zone_torque_deriv,
2050  const std::valarray<Eigen::VectorXd> &above_frac_deriv,
2051  unsigned locked_zone_ind,
2052  double *result
2053  ) const
2054  {
2055  unsigned global_zone_ind = zone_ind,
2056  num_zones = number_zones();
2057  if(&body == &__body2) global_zone_ind += __body1.number_zones();
2058  double above_frac =
2060  bool zone_is_locked = body.zone(zone_ind).locked();
2061  for(
2062  unsigned deriv_zone_ind = 0;
2063  deriv_zone_ind<num_zones;
2064  ++deriv_zone_ind
2065  ) {
2066  double &dest = (entry == Dissipation::PERIAPSIS
2067  ? result[deriv_zone_ind - 1]
2068  : result[deriv_zone_ind]);
2069  if(zone_is_locked)
2070  dest = (above_frac_deriv[deriv_zone_ind][locked_zone_ind]
2071  *
2072  (zone_z_torque_above - zone_z_torque_below));
2073  else dest = 0;
2074  if(std::abs(static_cast<int>(deriv_zone_ind
2075  -
2076  global_zone_ind)) <= 1) {
2077  double torque_deriv =
2078  zone_torque_deriv[deriv_zone_ind + 1 - global_zone_ind][2];
2079  if(zone_is_locked) dest += (above_frac * zone_torque_deriv[3][2]
2080  +
2081  (1.0 - above_frac) * torque_deriv);
2082  else dest += torque_deriv;
2083  }
2084  }
2085  }
2086 
2087  void BinarySystem::binary_jacobian(double *param_derivs,
2088  double *age_derivs) const
2089  {
2090  unsigned body1_zones = __body1.number_zones(),
2091  num_zones = body1_zones + __body2.number_zones(),
2092  num_locked_zones = (__body1.number_locked_zones()
2093  +
2095  num_param = 1 + 6 * num_zones - num_locked_zones;
2096  double dangmom_da = __orbital_angmom / (2.0 * __semimajor),
2097  dangmom_de = (-__eccentricity
2098  /
2099  (1.0 - std::pow(__eccentricity, 2))
2100  *
2102  std::valarray<double> orbit_power_deriv(num_param + 1),
2103  orbit_angmom_gain_deriv(num_param + 1);
2104  fill_orbit_power_deriv(orbit_power_deriv);
2105  fill_orbit_angmom_gain_deriv(orbit_angmom_gain_deriv);
2106  semimajor_jacobian(orbit_power_deriv,
2107  num_locked_zones == 0,
2108  param_derivs,
2109  age_derivs[0]);
2110  eccentricity_jacobian(orbit_power_deriv,
2111  orbit_angmom_gain_deriv,
2112  num_locked_zones == 0,
2113  param_derivs + num_param,
2114  age_derivs[1]);
2115  unsigned locked_zone_ind = 0;
2116  DissipatingBody *body = &__body1;
2117  unsigned body_zone_ind = 0;
2118  static Dissipation::QuantityEntry zone_deriv_list[]={
2122  };
2123  std::valarray<double> reference_periapsis_param_deriv(num_param);
2124  double reference_periapsis_age_deriv;
2125  const std::valarray<Eigen::VectorXd> *above_frac_deriv[] = {
2129  };
2130  for(unsigned zone_ind = 0; zone_ind < num_zones; ++zone_ind) {
2131  DissipatingZone &zone = body->zone(zone_ind);
2132  double sin_inc = std::sin(zone.inclination()),
2133  cos_inc = std::cos(zone.inclination());
2134  double &periapsis_age_deriv = (
2135  zone_ind==0
2136  ? reference_periapsis_age_deriv
2137  : age_derivs[1 + num_zones + zone_ind]
2138  );
2139  angle_evolution_age_deriv(*body,
2140  zone_ind,
2141  sin_inc,
2142  cos_inc,
2143  locked_zone_ind,
2144  age_derivs[2 + zone_ind],
2145  periapsis_age_deriv);
2146 
2147  double *periapsis_row = (
2148  zone_ind==0
2149  ? &reference_periapsis_param_deriv[0]
2150  : param_derivs+(1 + zone_ind+num_zones) * num_param
2151  );
2152 
2153  angle_evolution_orbit_deriv(Dissipation::SEMIMAJOR,
2154  dangmom_da,
2155  *body,
2156  body_zone_ind,
2157  sin_inc,
2158  cos_inc,
2159  locked_zone_ind,
2160  param_derivs[(2 + zone_ind) * num_param],
2161  periapsis_row[0]);
2162  angle_evolution_orbit_deriv(
2164  dangmom_de,
2165  *body,
2166  body_zone_ind,
2167  sin_inc, cos_inc,
2168  locked_zone_ind,
2169  param_derivs[(2 + zone_ind) * num_param + 1],
2170  periapsis_row[1]
2171  );
2172 
2173  Eigen::Vector3d zone_torque_above = body->nontidal_torque(zone_ind),
2174  zone_torque_below = zone_torque_above,
2175  orbit_torque = body->tidal_orbit_torque(zone);
2176  zone_torque_above += body->tidal_torque(zone_ind, true);
2177  zone_torque_below += body->tidal_torque(zone_ind, false);
2178  std::valarray<Eigen::Vector3d>
2179  zone_torque_deriv(zone.locked() ? 4 : 3),
2180  orbit_torque_deriv(num_zones);
2181  unsigned param_offset = (2 + zone_ind) * num_param + 2;
2182  for(unsigned deriv_ind = 0; deriv_ind < 3; ++deriv_ind) {
2183  Dissipation::QuantityEntry deriv = zone_deriv_list[deriv_ind];
2184  fill_zone_torque_deriv(deriv,
2185  *body,
2186  zone_ind,
2187  zone_torque_deriv);
2188  fill_orbit_torque_deriv(deriv,
2189  *body,zone_ind,
2190  orbit_torque_deriv);
2191  inclination_evolution_zone_derivs(deriv,
2192  *body, zone_ind,
2193  zone_torque_above[0],
2194  zone_torque_below[0],
2195  zone_torque_deriv,
2196  orbit_torque,
2197  orbit_torque_deriv,
2198  *(above_frac_deriv[deriv_ind]),
2199  sin_inc,
2200  cos_inc,
2201  locked_zone_ind,
2202  param_derivs + param_offset);
2203  periapsis_evolution_zone_derivs(deriv,
2204  *body,
2205  zone_ind,
2206  zone_torque_above[1],
2207  zone_torque_below[1],
2208  zone_torque_deriv,
2209  orbit_torque[1],
2210  orbit_torque_deriv,
2211  *(above_frac_deriv[deriv_ind]),
2212  sin_inc,
2213  cos_inc,
2214  locked_zone_ind,
2215  periapsis_row + 2);
2216  spin_angmom_evolution_zone_derivs(
2217  deriv,
2218  *body,
2219  zone_ind,
2220  zone_torque_above[2],
2221  zone_torque_below[2],
2222  zone_torque_deriv,
2223  *(above_frac_deriv[deriv_ind]),
2224  locked_zone_ind,
2225  param_derivs + param_offset + 2 * num_zones - 1
2226  );
2227  param_offset += (deriv == Dissipation::PERIAPSIS
2228  ? num_zones - 1
2229  : num_zones);
2230  }
2231 
2232  if(zone.locked()) ++locked_zone_ind;
2233  ++body_zone_ind;
2234  if(body_zone_ind == body->number_zones()) {
2235  body_zone_ind = 0;
2236  body = &__body2;
2237  }
2238  if(zone_ind != 0) {
2239  periapsis_age_deriv -= reference_periapsis_age_deriv;
2240  for(unsigned i = 0; i < num_param; ++i)
2241  periapsis_row[i] -= reference_periapsis_param_deriv[i];
2242  }
2243  }
2244 
2245  }
2246 #endif
2247 
2249  std::valarray<double> &orbit
2250  ) const
2251  {
2252  assert(__evolution_mode == Core::LOCKED_SURFACE_SPIN);
2253  assert(std::isnan(__semimajor));
2254  assert(std::isnan(__eccentricity));
2255 
2256  orbit.resize(__body1.number_zones() - 1);
2257  for(unsigned i = 1; i < __body1.number_zones(); ++i) {
2258  orbit[i - 1] = __body1.zone(i).angular_momentum();
2259  assert(__body1.zone(i).inclination() == 0);
2260  assert(__body1.zone(i).periapsis() == 0);
2261  }
2262  }
2263 
2264  void BinarySystem::fill_binary_orbit(std::valarray<double> &orbit) const
2265  {
2266  assert(__evolution_mode == Core::BINARY);
2267  assert(!std::isnan(__semimajor));
2268  assert(!std::isnan(__eccentricity));
2269  assert(__body1.zone(0).periapsis() == 0);
2270 
2271  orbit.resize(1 + 3 * number_zones() - number_locked_zones());
2272  if(number_locked_zones() == 0) orbit[0] = std::pow(__semimajor, 6.5);
2273  else orbit[0] = __semimajor;
2274  orbit[1] = __eccentricity;
2275  unsigned inclination_ind = 2,
2276  periapsis_ind = 2 + number_zones(),
2277  angmom_ind = periapsis_ind + number_zones() - 1;
2278  for(short body_ind = 0; body_ind < 2; ++body_ind) {
2279  DissipatingBody &body = (body_ind == 0 ? __body1 : __body2);
2280  for(
2281  unsigned zone_ind = 0;
2282  zone_ind < body.number_zones();
2283  ++zone_ind
2284  ) {
2285  DissipatingZone &zone = body.zone(zone_ind);
2286  orbit[inclination_ind++] = zone.inclination();
2287  if(body_ind || zone_ind)
2288  orbit[periapsis_ind++] = zone.periapsis();
2289  if(!zone.locked())
2290  orbit[angmom_ind++] = zone.angular_momentum();
2291  }
2292  }
2293  }
2294 
2295  void BinarySystem::fill_single_orbit(std::valarray<double> &orbit) const
2296  {
2297  assert(__evolution_mode == Core::SINGLE);
2298  assert(std::isnan(__semimajor));
2299  assert(std::isnan(__eccentricity));
2300  assert(__body1.zone(0).periapsis() == 0);
2301  assert(__body1.zone(0).inclination() == 0);
2302  assert(__body1.number_locked_zones() == 0);
2303 
2304  orbit.resize(3 * __body1.number_zones() - 2);
2305  unsigned inclination_ind = 0,
2306  periapsis_ind = __body1.number_zones() - 1,
2307  angmom_ind = periapsis_ind + __body1.number_zones() - 1;
2308  for(
2309  unsigned zone_ind = 0;
2310  zone_ind < __body1.number_zones();
2311  ++zone_ind
2312  ) {
2313  DissipatingZone &zone = __body1.zone(zone_ind);
2314  if(zone_ind) orbit[inclination_ind++] = zone.inclination();
2315  if(zone_ind) orbit[periapsis_ind++] = zone.periapsis();
2316  assert(!zone.locked());
2317  orbit[angmom_ind++] = zone.angular_momentum();
2318  }
2319  }
2320 
2321  int BinarySystem::configure(bool initialize,
2322  double age,
2323  double semimajor,
2324  double eccentricity,
2325  const double *spin_angmom,
2326  const double *inclination,
2327  const double *periapsis,
2329  {
2330 #ifndef NDEBUG
2331  if(initialize)
2332  std::cerr << "Initializing BinarySystem." << std::endl;
2333  if(evolution_mode != Core::BINARY) {
2334  assert(std::isnan(semimajor));
2335  assert(std::isnan(eccentricity));
2336  }
2337  if(evolution_mode == Core::LOCKED_SURFACE_SPIN) {
2338  assert(inclination == NULL);
2339  assert(periapsis == NULL);
2340  }
2341 
2342  std::cerr << "Configuring binary with a = " << semimajor
2343  << ", e = " << eccentricity
2344  << " in " << evolution_mode << " mode"
2345  << std::endl;
2346 #endif
2347 
2348  if(
2349  evolution_mode == Core::BINARY
2350  &&
2351  !(
2352  semimajor > 0.0
2353  &&
2354  eccentricity >= 0.0
2355  &&
2356  eccentricity < 1.0
2357  )
2358  )
2359  return GSL_EDOM;
2360 
2362  double m1 = __body1.mass(),
2363  m2 = __body2.mass();
2364  __age = age;
2367  __orbital_energy = Core::orbital_energy(m1, m2, semimajor);
2368  __orbital_angmom = Core::orbital_angular_momentum(m1,
2369  m2,
2370  semimajor,
2371  eccentricity);
2372 #ifndef NDEBUG
2373  std::cerr << "Configuring primary." << std::endl;
2374 #endif
2375  __body1.configure(initialize,
2376  age,
2377  m2,
2378  semimajor,
2379  eccentricity,
2380  spin_angmom,
2381  inclination,
2382  periapsis,
2383  evolution_mode == Core::LOCKED_SURFACE_SPIN,
2384  evolution_mode != Core::BINARY,
2385  true);
2386 
2387 
2388  if(evolution_mode == Core::BINARY) {
2389  unsigned offset = __body1.number_zones();
2390 #ifndef NDEBUG
2391  std::cerr << "Configuring secondary." << std::endl;
2392 #endif
2393  __body2.configure(initialize,
2394  age,
2395  m1,
2396  semimajor,
2397  eccentricity,
2398  spin_angmom + offset - __body1.number_locked_zones(),
2399  inclination + offset,
2400  periapsis + offset - 1);
2403  } else
2405 
2407 #ifndef NDEBUG
2408  // if(evolution_mode == BINARY)
2409  // assert(__semimajor > minimum_semimajor());
2410 #endif
2411  return 0;
2412  }
2413 
2414  int BinarySystem::configure(bool initialize,
2415  double age,
2416  const double *parameters,
2418  {
2420  double semimajor, eccentricity;
2421  const double *spin_angmom, *inclination, *periapsis;
2422  unsigned num_zones = number_zones();
2423  if(evolution_mode == Core::BINARY) {
2424  if(parameters[0] < 0) {
2425 #ifndef NDEBUG
2426  std::cerr << "At t = " << age << " param: ";
2427  for(
2428  unsigned i = 0;
2429  i < 3 * num_zones + 1 - number_locked_zones();
2430  ++i
2431  ) {
2432  if(i) std::cerr << ", ";
2433  std::cerr << parameters[i];
2434  }
2435  std::cerr << std::endl;
2436 #endif
2437  return GSL_EDOM;
2438  }
2440  semimajor = parameters[0];
2441  } else
2442  semimajor = std::pow(parameters[0], 1.0 / 6.5);
2443  eccentricity = parameters[1];
2444  inclination = parameters+2;
2445  } else {
2446  semimajor = eccentricity=Core::NaN;
2447  if(evolution_mode == Core::SINGLE) inclination = parameters;
2448  else inclination = NULL;
2449  }
2450 
2451  if(evolution_mode == Core::LOCKED_SURFACE_SPIN) {
2452  periapsis = NULL;
2453  spin_angmom = parameters;
2454  } else {
2455  assert(inclination != NULL);
2456 
2457  periapsis = inclination+num_zones;
2458  if(evolution_mode == Core::SINGLE) --periapsis;
2459  spin_angmom = periapsis + num_zones - 1;
2460  }
2461  return configure(initialize,
2462  age,
2463  semimajor,
2464  eccentricity,
2465  spin_angmom,
2466  inclination,
2467  periapsis,
2468  evolution_mode);
2469  }
2470 
2472  std::valarray<double> &orbit
2473  ) const
2474  {
2475  if(__evolution_mode == Core::LOCKED_SURFACE_SPIN)
2477  else if(__evolution_mode == Core::BINARY) fill_binary_orbit(orbit);
2478  else fill_single_orbit(orbit);
2479  return __evolution_mode;
2480  }
2481 
2482  double BinarySystem::above_lock_fraction(unsigned locked_zone_index,
2483  Dissipation::QuantityEntry entry,
2484  unsigned deriv_zone_index,
2485  bool secondary_radius)
2486  {
2487  if(zone_specific(entry)) {
2488  assert(entry == Dissipation::INCLINATION
2489  ||
2490  entry == Dissipation::PERIAPSIS
2491  ||
2492  entry == Dissipation::SPIN_ANGMOM
2493  ||
2495 
2496  if(entry == Dissipation::INCLINATION)
2497  return __above_lock_fractions_inclination_deriv[deriv_zone_index]
2498  [locked_zone_index];
2499  else if(entry == Dissipation::PERIAPSIS)
2500  return __above_lock_fractions_periapsis_deriv[deriv_zone_index]
2501  [locked_zone_index];
2502  else if(entry == Dissipation::SPIN_ANGMOM)
2503  return __above_lock_fractions_angmom_deriv[deriv_zone_index]
2504  [locked_zone_index];
2505  else return __above_lock_fractions_inertia_deriv[deriv_zone_index]
2506  [locked_zone_index];
2507  } else {
2508  if(entry == Dissipation::RADIUS && secondary_radius)
2509  return
2510  __above_lock_fractions_body2_radius_deriv[locked_zone_index];
2511  else return __above_lock_fractions[entry][locked_zone_index];
2512  }
2513  }
2514 
2516  const double *parameters,
2518  double *differential_equations)
2519  {
2520 #ifndef NDEBUG
2521  std::cerr << "Finding differential equations at t = " << age
2522  << " in " << evolution_mode
2523  << " mode, with orbit[0] = " << parameters[0]
2524  << std::endl;
2525 #endif
2526  int status = configure(false, age, parameters, evolution_mode);
2527  if(status != GSL_SUCCESS) return status;
2528  __semimajor_rate = __eccentricity_rate = Core::NaN;
2529  switch(evolution_mode) {
2530  case Core::LOCKED_SURFACE_SPIN :
2532  differential_equations
2533  );
2534  case Core::SINGLE :
2536  differential_equations
2537  );
2538  case Core::BINARY :
2539  return binary_differential_equations(differential_equations);
2540  default :
2542  "Evolution mode other than LOCKED_SURFACE_SPIN, SINGLE or "
2543  "BINARY encountered in "
2544  "BinarySystem::differential_equations!"
2545  );
2546  }
2547  }
2548 
2549 #ifdef ENABLE_DERIVATIVES
2550  int BinarySystem::jacobian(double age,
2551  const double *parameters,
2553  double *param_derivs,
2554  double *age_derivs)
2555  {
2556  configure(false, age, parameters, evolution_mode);
2557  switch(evolution_mode) {
2558  case Core::LOCKED_SURFACE_SPIN :
2559  locked_surface_jacobian(param_derivs, age_derivs);
2560  return 0;
2561  case Core::SINGLE : single_body_jacobian(param_derivs, age_derivs);
2562  return 0;
2563  case Core::BINARY : binary_jacobian(param_derivs, age_derivs);
2564  return 0;
2565  default : throw Core::Error::BadFunctionArguments(
2566  "Evolution mode other than LOCKED_SURFACE_SPIN, "
2567  "SINGLE or BINARY encountered in "
2568  "BinarySystem::jacobian!"
2569  );
2570  }
2571  }
2572 #endif
2573 
2574  void BinarySystem::initialize_locks(double sync_precision)
2575  {
2576  lock_scenario_type lock_candidates = find_synchronized_zones(
2577  sync_precision
2578  );
2579  if(lock_candidates.empty())
2580  return;
2581 
2582  lock_scenario_type lock_scenario;
2583  if(!explore_lock_scenarios(lock_candidates.begin(),
2584  lock_candidates.size(),
2585  lock_scenario))
2586  throw Core::Error::Runtime("Failed to find viable lock scenario!");
2587 #ifndef NDEBUG
2588  std::cerr << "Setting lock scenario: " << std::endl;
2590  std::cerr,
2592  std::vector<bool>(__selected_lock_scenario.size(), true),
2593  true
2594  );
2596 #endif
2597  }
2598 
2599  void BinarySystem::check_for_lock(int orbital_freq_mult,
2600  int spin_freq_mult,
2601  unsigned short body_index,
2602  unsigned zone_index,
2603  short direction)
2604  {
2605  DissipatingBody &body = (body_index ? __body2 : __body1);
2606 
2607  assert(body_index <= 1);
2608  assert(zone_index < body.number_zones());
2609  assert(spin_freq_mult);
2610 
2611  double original_angmom = body.zone(zone_index).angular_momentum();
2612  body.lock_zone_spin(zone_index, orbital_freq_mult, spin_freq_mult);
2613  unsigned num_zones = number_zones(),
2614  num_locked_zones = number_locked_zones(),
2615  locked_zone_ind = 0;
2616  std::vector<double> spin_angmom(num_zones - num_locked_zones),
2617  inclinations(num_zones),
2618  periapses(num_zones);
2619  for(unsigned zone_ind = 0; zone_ind < num_zones; ++zone_ind) {
2620  DissipatingZone &zone=(zone_ind < __body1.number_zones()
2621  ? __body1.zone(zone_ind)
2622  : __body2.zone(zone_ind
2623  -
2624  __body1.number_zones()));
2625  if(zone.locked()) ++locked_zone_ind;
2626  else spin_angmom[zone_ind-locked_zone_ind] = (
2627  zone.angular_momentum()
2628  );
2629  inclinations[zone_ind] = zone.inclination();
2630  if(zone_ind) periapses[zone_ind - 1] = zone.periapsis();
2631  }
2632  assert(locked_zone_ind == num_locked_zones);
2633  configure(false,
2634  __age,
2635  __semimajor,
2637  &(spin_angmom[0]),
2638  &(inclinations[0]),
2639  &(periapses[0]),
2640  Core::BINARY);
2641  DissipatingZone &locked_zone = body.zone(zone_index);
2644  [locked_zone.locked_zone_index()];
2645 #ifndef NDEBUG
2646  std::cerr << "Holding lock requires above lock fraction of: "
2647  << above_lock_fraction
2648  << std::endl;
2649 #endif
2650  if(above_lock_fraction > 0 && above_lock_fraction < 1) return;
2651  std::vector<double>::iterator check_zone_dest = (
2652  spin_angmom.begin()
2653  +
2654  (
2655  zone_index
2656  +
2657  (body_index ? __body1.number_zones() : 0)
2658  -
2659  locked_zone.locked_zone_index()
2660  )
2661  );
2662  spin_angmom.insert(check_zone_dest, original_angmom);
2663  if(direction == 0)
2664  throw Core::Error::Runtime(
2665  "Crossing spin-orbit synchronization with unspecified direction!"
2666  );
2667  body.unlock_zone_spin(zone_index, direction);
2668  configure(false,
2669  __age,
2670  __semimajor,
2672  &(spin_angmom[0]),
2673  &(inclinations[0]),
2674  &(periapses[0]),
2675  Core::BINARY);
2676  }
2677 
2678  double BinarySystem::minimum_separation(bool deriv) const
2679  {
2680  double rroche = (
2681  __body2.radius()
2682  ?(
2683  2.44
2684  *
2685  __body2.radius()
2686  *
2687  std::pow(__body1.mass() / __body2.mass(), 1.0 / 3.0)
2688  )
2689  :(
2690  0.0
2691  )
2692  );
2693  if(rroche > __body1.radius()) {
2694  if(deriv) return rroche * __body2.radius(1) / __body2.radius();
2695  else return rroche;
2696  } else {
2697  return __body1.radius(deriv ? 1 : 0);
2698  }
2699  }
2700 
2702  {
2703 #ifndef NDEBUG
2704  std::cerr << "Handling secondary death!" << std::endl;
2705 #endif
2706  unsigned num_zones = __body1.number_zones();
2707  std::valarray<double> spin_angmom(num_zones),
2708  inclination(num_zones - 1),
2709  periapsis(num_zones - 1);
2710  DissipatingZone &old_surface_zone = __body1.zone(0);
2711  double old_surface_inclination = old_surface_zone.inclination(),
2712  sin_inc = std::sin(old_surface_inclination),
2713  cos_inc = std::cos(old_surface_inclination),
2714  old_surface_angmom = old_surface_zone.angular_momentum(),
2715  angmom = std::sqrt(
2716  std::pow(old_surface_angmom + __orbital_angmom * cos_inc, 2)
2717  +
2718  std::pow(__orbital_angmom * sin_inc, 2)
2719  ),
2720  new_surface_inclination=std::atan2(__orbital_angmom*sin_inc,
2721  old_surface_angmom
2722  +
2723  __orbital_angmom*cos_inc);
2724 
2725  spin_angmom[0] = angmom;
2726  // assert(num_zones == 2);
2727  if(old_surface_zone.locked()) __body1.unlock_zone_spin(0, 1);
2728  for(unsigned zone_ind = 1; zone_ind < num_zones; ++zone_ind) {
2729  DissipatingZone &zone = __body1.zone(zone_ind);
2730  assert(zone.periapsis() == 0);
2731  inclination[zone_ind - 1]=std::abs(zone.inclination()
2732  -
2733  old_surface_inclination
2734  +
2735  new_surface_inclination);
2736  periapsis[zone_ind - 1] = 0;
2737  spin_angmom[zone_ind] = zone.angular_momentum();
2738  if(zone.locked()) __body1.unlock_zone_spin(zone_ind, 1);
2739  }
2740  configure(false,
2741  __age,
2742  Core::NaN,
2743  Core::NaN,
2744  &spin_angmom[0],
2745  &inclination[0],
2746  &periapsis[0],
2747  Core::SINGLE);
2748  __body1.spin_jumped();
2749  }
2750 
2751  void BinarySystem::release_lock(unsigned locked_zone_index, short direction)
2752  {
2753  DissipatingBody *body;
2754  if(locked_zone_index < __body1.number_locked_zones()) {
2755  body = &__body1;
2756  for(
2757  unsigned zone_ind = 0;
2758  zone_ind < __body2.number_locked_zones();
2759  ++zone_ind
2760  )
2761  if(__body2.zone(zone_ind).locked())
2762  --__body2.zone(zone_ind).locked_zone_index();
2763 
2764  } else {
2765  body = &__body2;
2766  locked_zone_index -= __body1.number_locked_zones();
2767  }
2768  unsigned zone_ind = 0;
2769  while(true){
2770  if(body->zone(zone_ind).locked()) {
2771  if(locked_zone_index == 0) break;
2772  else --locked_zone_index;
2773  }
2774  assert(zone_ind < body->number_zones());
2775  ++zone_ind;
2776  }
2777  body->unlock_zone_spin(zone_ind, direction);
2778  for(; zone_ind < body->number_zones(); ++zone_ind)
2779  if(body->zone(zone_ind).locked())
2780  --body->zone(zone_ind).locked_zone_index();
2781  }
2782 
2784  {
2791  }
2792 
2794  {
2795  __semimajor_evolution.clear();
2796  __eccentricity_evolution.clear();
2801  }
2802 
2803  void BinarySystem::rewind_evolution(unsigned nsteps)
2804  {
2805  for(unsigned i = 0; i < nsteps; ++i) {
2806  __semimajor_evolution.pop_back();
2807  __eccentricity_evolution.pop_back();
2808  __semimajor_rate_evolution.pop_back();
2809  __eccentricity_rate_evolution.pop_back();
2810  }
2811  __body1.rewind_evolution(nsteps);
2812  __body2.rewind_evolution(nsteps);
2813  }
2814 
2816  {
2818  if(__evolution_mode == Core::BINARY)
2819  (*result) |= new SecondaryDeathCondition(*this);
2820  (*result) |= __body1.stopping_conditions(*this, true);
2821  if(__evolution_mode == Core::BINARY)
2822  (*result) |= __body2.stopping_conditions(*this, false);
2823  return result;
2824  }
2825 
2827  {
2830  }
2831 
2833  {
2834  return std::min(__body1.next_stop_age(),
2836  }
2837 
2839  {
2840 #ifndef NDEBUG
2841  int result = -1;
2842 #endif
2843  DissipatingBody *body = &__body1;
2844  while(true) {
2845  for(
2846  unsigned zone_ind = 0;
2847  zone_ind < body->number_zones();
2848  ++zone_ind
2849  ) {
2850  DissipatingZone &zone = body->zone(zone_ind);
2851  if(zone.dissipative()) {
2852 #ifndef NDEBUG
2853  if(result < 1)
2854  result = zone.expansion_order();
2855  else
2856  assert(static_cast<unsigned>(result)
2857  ==
2858  zone.expansion_order());
2859 #else
2860  return zone.expansion_order();
2861 #endif
2862  }
2863  }
2864  if(body == &__body2) break;
2865  body = &__body2;
2866  }
2867 #ifndef NDEBUG
2868  if(result >= 1)
2869  return result;
2870 #endif
2872  "System contains no dissipative zones! "
2873  "Run as single objects instead!"
2874  );
2875  }
2876 
2877 } //End Evolve namespace.
virtual void add_to_evolution()
Appends the state defined by last configure(), to the evolution.
double __age
The present age of the stellar system in Gyrs.
Definition: BinarySystem.h:81
unsigned locked_zone_index() const
The index of this zone in the list of locked zones (valid only if locked).
void initialize_locks(double sync_precision)
Identify and lock all zones within precision of a lock that can hold the lock at the current configur...
virtual void rewind_evolution(unsigned nsteps)
Discards the last steps from the evolution.
virtual void reset_evolution()
Discards all evolution.
int single_body_differential_equations(double *evolution_rates) const
Differential equations for the rotation of the zones of body 0 if no other body is present...
Eigen::VectorXd __above_lock_fractions_body2_radius_deriv
The derivative of the above lock fractions w.r.t. to the radius of the secondardy.
Definition: BinarySystem.h:148
Function arguments do not satisfy some requirement.
Definition: Error.h:73
RADIUS
The derivative w.r.t. the radius of the body in .
const std::list< double > & eccentricity_evolution() const
The tabulated evolution of the eccentricity so far.
virtual void add_to_evolution()
Appends the state defined by last configure(), to the evolution.
Eigen::Vector3d __orbit_torque
The torque on the orbit in the coordinate system of the outermost zone of the first body...
Definition: BinarySystem.h:112
const SpinOrbitLockInfo & lock_monitored(bool other=false) const
Return either of the locks being monitored.
void find_locked_zones()
Fills the __locked_zones list.
double angular_momentum() const
The angular momentum of the given zone in .
double __semimajor_rate
The current rate at which the semiamjor axis is changing.
Definition: BinarySystem.h:105
unsigned number_locked_zones() const
The number of zones currently in a spin-orbit lock.
Satisfied when some multiples of the orbit and stellar rotation are synchronized. ...
virtual unsigned expansion_order() const
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.
int binary_differential_equations(double *differential_equations) const
The differential equations for a system with both bodies present.
std::list< double > __eccentricity_evolution
The evolution of the eccentricity recorded by add_to_evolution() so far.
Definition: BinarySystem.h:65
void set_lock_scenario(const lock_scenario_type &lock_scenario)
Lock exactly the specified zones and configure the system.
SPIN_FREQUENCY
The derivative w.r.t. the spin frequency of a dissipating zone.
std::list< double > __semimajor_evolution
The evolution of the semimajor axis recorded by add_to_evolution() so far.
Definition: BinarySystem.h:65
A base class for any body contributing to tidal dissipation.
double tidal_power(unsigned zone_index, bool above, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV) const
Tidal power dissipated in the given zone.
virtual double minimum_separation(bool deriv=false) const
Smallest semimajor axis at which the secondary can survive for the latest system configuration.
int orbital_frequency_multiplier() const
The multiplier in front of the orbital frequency in the lock.
int locked_surface_differential_equations(double *evolution_rates) const
Differential equations for the rotation of the zones of body 0 with the topmost zone rotating with a ...
void fill_above_lock_fractions_deriv()
Fills the __above_lock_fractinos_*_deriv members.
double periapsis(bool evolution_rate=false) const
The argument of periapsis of this zone minus the reference zone&#39;s.
virtual const DissipatingZone & zone(unsigned zone_index) const =0
A modifiable reference to one of the body&#39;s zones.
double __semimajor
The current semimajor axis.
Definition: BinarySystem.h:81
void check_for_lock(int orbital_freq_mult, int spin_freq_mult, unsigned short body_index, unsigned zone_index, short direction)
Check if a spin-orbit lock can be held and updates the system as necessary to calculate subsequent ev...
std::valarray< Eigen::VectorXd > __above_lock_fractions
The above lock fractinos for the locked zones and their derivatives.
Definition: BinarySystem.h:125
SEMIMAJOR
The derivative w.r.t. the semimajor axis in AU.
void unlock_zone_spin(unsigned zone_index, short direction)
Releases the given zone from a spin-orbit lock.
MOMENT_OF_INERTIA
virtual void reached_critical_age(double)
Change the body as necessary at the given age.
virtual double next_stop_age() const
The next age when the evolution needs to be stopped for a system change.
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 above_lock_fraction(unsigned locked_zone_index, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, unsigned deriv_zone_index=0, bool secondary_radius=false)
The fraction of an infinitesimal timestep that a zone spends spinning faster than the lock it is in...
Any runtime error.
Definition: Error.h:61
Orientations of zones of bodies in a binary system.
const std::list< double > & semimajor_evolution() const
The tabulated evolution of the semimajor axis so far.
lock_scenario_type __selected_lock_scenario
The lock/unlock configuration selected by the last call to.
Definition: BinarySystem.h:178
bool test_synchronized_unlocked_zone(unsigned test_zone_ind)
NUM_DERIVATIVES
The total number of derivatives supported.
double __orbit_power
The rate at which the orbit gains energy due to tides.
Definition: BinarySystem.h:81
void fill_orbit_angmom_gain_deriv(std::valarray< double > &orbit_angmom_gain_deriv) const
Computes the derivatives w.r.t. the evolution quantities of the orbit angular momentum gain...
int differential_equations(double age, const double *parameters, Core::EvolModeType evolution_mode, double *differential_equations)
The differential equation and jacobian for the evolution of the system.
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...
virtual void reset_evolution()
Resets the evolution of the system.
double mass() const
The mass of the body (constant with age).
void add_body_rate_deriv(const DissipatingBody &body, VALUE_TYPE(DissipatingBody::*func)(Dissipation::QuantityEntry, unsigned, const Eigen::VectorXd &) const, std::valarray< VALUE_TYPE > &orbit_rate_deriv, unsigned offset) const
Adds the derivatives of a rate by which the orbit is changing due to tidal interactions with a single...
double age() const
Returns the present age of the system in Gyr.
Definition: BinarySystem.h:914
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...
void fill_locked_surface_orbit(std::valarray< double > &orbit) const
Implements fill_orbit() for LOCKED_SURFACE_SPIN evolution mode.
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()).
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...
int spin_frequency_multiplier() const
The multiplier in front of the spin frequency in the lock.
virtual unsigned expansion_order() const
std::valarray< Eigen::VectorXd > __above_lock_fractions_angmom_deriv
The derivatives of the above lock fractions w.r.t. the angular momenta of the zones.
Definition: BinarySystem.h:125
unsigned number_locked_zones() const
How many zones on either body are currently locked.
Definition: BinarySystem.h:929
std::list< double > __semimajor_rate_evolution
The rate at which the semimajor axis evolved at each.
Definition: BinarySystem.h:65
DissipatingBody & __body2
The second body in the system.
Definition: BinarySystem.h:163
virtual void reached_critical_age(double age)
Change the system as necessary at the given age.
void calculate_above_lock_fractions(Eigen::VectorXd &fractions, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, bool body1_deriv=true)
Calculates the fraction of a timestep above spin-orbit lock for all locked zones. ...
double __eccentricity
The current eccentricity.
Definition: BinarySystem.h:81
virtual double next_stop_age() const
The next age when the evolution needs to be stopped for a change in one of the bodies.
void fill_binary_orbit(std::valarray< double > &orbit) const
Implements fill_orbit() for BINARY evolution mode.
Satisfied when the planet enters below either the roche sphere or the stellar photosphere.
std::list< unsigned > __locked_zones
A list of indices of locked zones.
Definition: BinarySystem.h:118
double __orbit_angmom_gain
The rate at which the orbit gains angular momentum due to tides.
Definition: BinarySystem.h:81
void fill_single_orbit(std::valarray< double > &orbit) const
Implements fill_orbit() for SINGLE evolution mode.
std::list< double > __eccentricity_rate_evolution
The rate at which the eccentricity evolved at each.
Definition: BinarySystem.h:65
Core::EvolModeType __evolution_mode
The evolution mode from the last call to configure();.
Definition: BinarySystem.h:115
virtual void secondary_died()
Update the system to account for the death of the secondary.
Eigen::ColPivHouseholderQR< Eigen::MatrixXd > __above_lock_fractions_decomp
The matrix that defines the problem to solve for __above_lock_fractions.
Definition: BinarySystem.h:155
bool test_lock_scenario(const lock_scenario_type &lock_scenario, bool first_scenario)
Return true iff all locks in the given scenario will be maintained and synchronized unlocked zones wi...
double periapsis_evolution(const Eigen::Vector3d &orbit_torque, const Eigen::Vector3d &zone_torque, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, const Eigen::Vector3d &orbit_torque_deriv=Eigen::Vector3d(), const Eigen::Vector3d &zone_torque_deriv=Eigen::Vector3d())
The rate at which the periapsis of the orbit/reference zone in this zone&#39;s equatorial plane is changi...
ECCENTRICITY
The derivative w.r.t. the eccentricity.
double orbital_frequency(bool semimajor_deriv=false) const
The current orbital frequency [Rad/day] (calculated from the semimajor axis).
Definition: BinarySystem.h:940
void set_evolution_rates(double angular_momentum, double inclination, double periapsis)
Set evolution rates for this zone&#39;s properties for storing in the eveloution.
void set_reference_zone_angmom(double reference_angmom)
Defines the angular momentum of the reference zone for single body evolution.
std::valarray< Eigen::VectorXd > __above_lock_fractions_inclination_deriv
The derivatives of the above lock fractions w.r.t. the inclinations of the zones. ...
Definition: BinarySystem.h:125
AGE
The derivative w.r.t. age, excluding the dependence through the body&#39;s radius and the moments of iner...
void fill_orbit_power_deriv(std::valarray< double > &orbit_power_deriv) const
Computes the derivatives w.r.t. the evolution quantities of the orbit energy gain.
Defines the BinarySystem class.
virtual void rewind_evolution(unsigned nsteps)
Discards the last steps from the evolution.
lock_scenario_type find_synchronized_zones(double precision)
Return a list of the zone indices that are spin-orbit synchronized with some tidal term (locked or un...
Core::EvolModeType evolution_mode()
The evolution mode of last call to configure().
void update_above_lock_fractions()
Solves for and sets the above lock fractions and their derivatives.
virtual bool locked(int orbital_frequency_multiplier, int spin_frequency_multiplier) const
Should return true iff the given term is presently locked.
double inclination(bool evolution_rate=false) const
The angle between the angular momenta of the zone and the orbit.
NO_DERIV
The quantity itself, undifferentiated.
EvolModeType
The various evolution modes.
Definition: Common.h:42
ORBITAL_FREQUENCY
The derivative w.r.t. the orbital frequency.
bool explore_lock_scenarios(lock_scenario_type::const_iterator next_synchronized_zone, unsigned num_synchronized_zones, lock_scenario_type &lock_scenario, bool first_scenario=true)
Iterate over all possible lock scenarios until a valid one is found.
std::valarray< Eigen::VectorXd > __above_lock_fractions_inertia_deriv
The derivatives of the above lock fractions w.r.t. the moments of inertia of the zones.
Definition: BinarySystem.h:125
void lock_zone_spin(unsigned zone_index, int orbital_frequency_multiplier, int spin_frequency_multiplier)
double semimajor() const
The current semimajor axis of the system.
Definition: BinarySystem.h:933
SPIN_ANGMOM
The derivative w.r.t. the spin angular momentum in .
double inclination_evolution(const Eigen::Vector3d &orbit_torque, const Eigen::Vector3d &zone_torque, Dissipation::QuantityEntry entry=Dissipation::NO_DERIV, const Eigen::Vector3d &orbit_torque_deriv=Eigen::Vector3d(), const Eigen::Vector3d &zone_torque_deriv=Eigen::Vector3d())
The rate at which the inclination between this zone and the orbit is changing.
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()).
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()).
DissipatingBody & __body1
The first body in the system.
Definition: BinarySystem.h:163
Eigen::VectorXd above_lock_fractions_deriv(Dissipation::QuantityEntry entry, DissipatingBody &body, unsigned zone_index)
Calculates derivatives of the above lock fractions w.r.t. quantities of non-surface zones...
std::valarray< Eigen::VectorXd > __above_lock_fractions_periapsis_deriv
The derivatives of the above lock fractions w.r.t. the periapses of the zones.
Definition: BinarySystem.h:125
A class combining the the outputs of multiple stopping conditions.
double __eccentricity_rate
The current rate at which the eccentricity is changing.
Definition: BinarySystem.h:105
virtual bool dissipative() const =0
Should return true iff the zone has some non-zero dissipation.
unsigned number_zones() const
The total number of zones in both system bodies.
Definition: BinarySystem.h:923
void above_lock_problem_deriv_correction(Dissipation::QuantityEntry entry, bool body1_deriv, Eigen::MatrixXd &matrix, Eigen::VectorXd &rhs) const
Makes corrections to the matrix and RHS to accomodate the given derivative for the linear problem tha...
virtual bool can_lock() const =0
Should return true iff the zone&#39;s dissipation is discontinuous at zero frequency. ...
virtual void release_lock(unsigned locked_zone_index, short direction)
Releases the lock to one of the locked zones.
void unlock_all_zones(const std::vector< short > &unlock_directions, const std::vector< double > &original_angmom=std::vector< double >())
Unlock all zones, restoring their original spins if known.
double __orbital_angmom
The current orbital angular momentum.
Definition: BinarySystem.h:81
double eccentricity() const
The current eccentricity of the system.
Definition: BinarySystem.h:936
virtual CombinedStoppingCondition * stopping_conditions()
Conditions detecting the next possible doscontinuity in the evolution.
virtual int configure(bool initialize, double age, double semimajor, double eccentricity, const double *spin_angmom, const double *inclination, const double *periapsis, Core::EvolModeType evolution_mode)
Sets the current state of the system.
double __orbital_energy
The current orbital energy.
Definition: BinarySystem.h:81
void fill_orbit_torque_and_power()
Update the values of ::__orbit_power, ::__orbit_torque, and ::__orbit_angmom_gain.
void describe_lock_scenario(std::ostream &os, const lock_scenario_type &lock_scenario, const std::vector< bool > passed, bool add_header=false)
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.
Core::EvolModeType fill_orbit(std::valarray< double > &orbit) const
Fills an array with the parameters expected by differential_equations() and jacobian(), returning the evolution mode.
double spin_frequency() const
The spin frequency of the given zone.
Defines a lock between the spin of a dissipating body and the orbit.
virtual unsigned number_zones() const =0
The number of zones the body consists of.