8 #ifndef __BINARY_SYSTEM_H 9 #define __BINARY_SYSTEM_H 11 #include "../Core/SharedLibraryExportMacros.h" 15 #include "../Core/AstronomicalConstants.h" 16 #include "../Core/Common.h" 17 #include "../Core/OrbitalExpressions.h" 18 #include "../Core/Error.h" 19 #include <gsl/gsl_errno.h> 20 #include <gsl/gsl_odeiv2.h> 21 #include <gsl/gsl_siman.h> 65 __semimajor_evolution,
69 __eccentricity_evolution,
77 __eccentricity_rate_evolution;
125 std::valarray<Eigen::VectorXd> __above_lock_fractions,
129 __above_lock_fractions_inclination_deriv,
137 __above_lock_fractions_inertia_deriv,
141 __above_lock_fractions_angmom_deriv;
154 Eigen::ColPivHouseholderQR<Eigen::MatrixXd>
172 std::list< std::tuple<unsigned, short, double> >
182 void find_locked_zones();
198 lock_scenario_type find_synchronized_zones(
208 void set_lock_scenario(
212 const lock_scenario_type &lock_scenario
216 void unlock_all_zones(
218 const std::vector<short> &unlock_directions,
224 const std::vector<double> &original_angmom = std::vector<double>()
229 bool test_synchronized_unlocked_zone(
232 unsigned test_zone_ind
236 void describe_lock_scenario(
242 const lock_scenario_type &lock_scenario,
246 const std::vector<bool> passed,
249 bool add_header=
false 255 bool test_lock_scenario(
258 const lock_scenario_type &lock_scenario
261 ,
bool first_scenario
269 bool explore_lock_scenarios(
272 lock_scenario_type::const_iterator next_synchronized_zone,
275 unsigned num_synchronized_zones,
278 lock_scenario_type &lock_scenario
282 ,
bool first_scenario =
true 291 int locked_surface_differential_equations(
293 double *evolution_rates
296 #ifdef ENABLE_DERIVATIVES 297 void locked_surface_jacobian(
304 double *param_derivs,
316 int single_body_differential_equations(
319 double *evolution_rates
322 #ifdef ENABLE_DERIVATIVES 323 void fill_single_body_jacobian(
327 double *inclination_param_derivs,
330 double *periapsis_param_derivs,
334 double *angmom_param_derivs,
338 double *inclination_age_derivs,
342 double *periapsis_age_derivs,
346 double *angmom_age_derivs
354 void single_body_jacobian(
356 double *param_derivs,
374 double semimajor_evolution(
383 double orbit_power_deriv=Core::NaN
396 double eccentricity_evolution(
403 double orbit_angmom_gain,
410 double orbit_power_deriv=Core::NaN,
415 double orbit_angmom_gain_deriv=Core::NaN,
419 bool semimajor_deriv=
true 424 void above_lock_problem_deriv_correction(
426 Dissipation::QuantityEntry entry,
433 Eigen::MatrixXd &matrix,
446 void calculate_above_lock_fractions(
448 Eigen::VectorXd &fractions,
458 bool body1_deriv=
true 463 Eigen::VectorXd above_lock_fractions_deriv(
466 Dissipation::QuantityEntry entry,
478 void fill_above_lock_fractions_deriv();
486 void update_above_lock_fractions();
492 void fill_orbit_torque_and_power();
499 int binary_differential_equations(
502 double *differential_equations
512 template<
typename VALUE_TYPE>
513 void add_body_rate_deriv(
521 const Eigen::VectorXd &)
const,
525 std::valarray<VALUE_TYPE> &orbit_rate_deriv,
534 void fill_orbit_power_deriv(
536 std::valarray<double> &orbit_power_deriv
541 void fill_orbit_angmom_gain_deriv(
543 std::valarray<double> &orbit_angmom_gain_deriv
546 #ifdef ENABLE_DERIVATIVES 547 void semimajor_jacobian(
551 const std::valarray<double> &orbit_power_deriv,
558 double *param_derivs,
566 void eccentricity_jacobian(
569 const std::valarray<double> &orbit_power_deriv,
573 const std::valarray<double> &orbit_angmom_gain_deriv,
580 double *param_derivs,
589 void angle_evolution_age_deriv(
605 unsigned locked_zone_ind,
618 void angle_evolution_orbit_deriv(
621 Dissipation::QuantityEntry entry,
640 unsigned locked_zone_ind,
650 void fill_orbit_torque_deriv(
654 Dissipation::QuantityEntry entry,
668 std::valarray<Eigen::Vector3d> &orbit_torque_deriv
673 void fill_zone_torque_deriv(
677 Dissipation::QuantityEntry entry,
692 std::valarray<Eigen::Vector3d> &zone_torque_deriv
697 void inclination_evolution_zone_derivs(
701 Dissipation::QuantityEntry entry,
711 double zone_x_torque_above,
715 double zone_x_torque_below,
718 const std::valarray<Eigen::Vector3d> &zone_torque_deriv,
722 const Eigen::Vector3d &orbit_torque,
725 const std::valarray<Eigen::Vector3d> &orbit_torque_deriv,
728 const std::valarray<Eigen::VectorXd> &above_frac_deriv,
738 unsigned locked_zone_ind,
746 void periapsis_evolution_zone_derivs(
748 Dissipation::QuantityEntry entry,
758 double zone_y_torque_above,
762 double zone_y_torque_below,
765 const std::valarray<Eigen::Vector3d> &zone_torque_deriv,
768 double orbit_y_torque,
771 const std::valarray<Eigen::Vector3d> &orbit_torque_deriv,
774 const std::valarray<Eigen::VectorXd> &above_frac_deriv,
784 unsigned locked_zone_ind,
790 void spin_angmom_evolution_zone_derivs(
792 Dissipation::QuantityEntry entry,
802 double zone_z_torque_above,
806 double zone_z_torque_below,
809 const std::valarray<Eigen::Vector3d> &zone_torque_deriv,
812 const std::valarray<Eigen::VectorXd> &above_frac_deriv,
816 unsigned locked_zone_ind,
823 void binary_jacobian(
825 double *param_derivs,
834 void fill_locked_surface_orbit(std::valarray<double> &orbit)
const;
837 void fill_binary_orbit(std::valarray<double> &orbit)
const;
840 void fill_single_orbit(std::valarray<double> &orbit)
const;
854 const std::string &system_name=
"" 863 const std::string
get_name()
const {
return __name;}
866 virtual int configure(
881 const double *spin_angmom,
886 const double *inclination,
891 const double *periapsis,
907 const double *parameters,
914 double age()
const {
return __age;}
924 {
return (__evolution_mode == Core::BINARY
925 ? __body1.number_zones() + __body2.number_zones()
926 : __body1.number_zones());}
930 {
return __body1.number_locked_zones() + __body2.number_locked_zones();}
942 return Core::orbital_angular_velocity(__body1.mass(),
954 std::valarray<double> &orbit
960 double above_lock_fraction(
962 unsigned locked_zone_index,
972 unsigned deriv_zone_index=0,
976 bool secondary_radius=
false 986 int differential_equations(
1016 const double *parameters,
1024 double *differential_equations
1027 #ifdef ENABLE_DERIVATIVES 1036 const double *parameters,
1043 double *param_derivs,
1053 void initialize_locks(
1057 double sync_precision
1062 void check_for_lock(
1064 int orbital_freq_mult,
1071 unsigned short body_index,
1074 unsigned zone_index,
1087 virtual double minimum_separation(
1099 virtual void secondary_died();
1102 virtual void release_lock(
1105 unsigned locked_zone_index,
1113 virtual void add_to_evolution();
1116 virtual void reset_evolution();
1119 virtual void rewind_evolution(
1134 virtual void reached_critical_age(
double age);
1138 virtual double next_stop_age()
const;
1142 {
return __semimajor_evolution;}
1146 {
return __semimajor_rate_evolution;}
1150 {
return __eccentricity_evolution;}
1154 {
return __eccentricity_rate_evolution;}
1159 unsigned new_expansion_order
1162 __body1.change_expansion_order(new_expansion_order, *
this,
true);
1163 __body2.change_expansion_order(new_expansion_order, *
this,
false);
1168 virtual unsigned expansion_order()
const;
Eigen::VectorXd __above_lock_fractions_body2_radius_deriv
The derivative of the above lock fractions w.r.t. to the radius of the secondardy.
const std::list< double > & eccentricity_evolution() const
The tabulated evolution of the eccentricity so far.
const std::string get_name() const
Returns the name of the system.
Eigen::Vector3d __orbit_torque
The torque on the orbit in the coordinate system of the outermost zone of the first body...
double __semimajor_rate
The current rate at which the semiamjor axis is changing.
const DissipatingBody & primary() const
Returns the primary body in the system (const).
A base class for any body contributing to tidal dissipation.
std::string __name
The name of the binary system (e.g. "HAT-P-20")
double __semimajor
The current semimajor axis.
const std::list< double > & eccentricity_evolution_rate() const
The tabulated evolution of the eccentricity so far.
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.
NUM_DERIVATIVES
The total number of derivatives supported.
double age() const
Returns the present age of the system in Gyr.
unsigned number_locked_zones() const
How many zones on either body are currently locked.
const DissipatingBody & secondary() const
Returns the secondary body in the system (const).
std::list< double > __semimajor_rate_evolution
The rate at which the semimajor axis evolved at each.
DissipatingBody & __body2
The second body in the system.
BinarySystem(DissipatingBody &body1, DissipatingBody &body2, const std::string &system_name="")
Construct a binary system.
std::list< unsigned > __locked_zones
A list of indices of locked zones.
Core::EvolModeType __evolution_mode
The evolution mode from the last call to configure();.
Eigen::ColPivHouseholderQR< Eigen::MatrixXd > __above_lock_fractions_decomp
The matrix that defines the problem to solve for __above_lock_fractions.
double orbital_frequency(bool semimajor_deriv=false) const
The current orbital frequency [Rad/day] (calculated from the semimajor axis).
Declares a stopping condition class monitoring for the death of the secondary object.
Core::EvolModeType evolution_mode()
The evolution mode of last call to configure().
NO_DERIV
The quantity itself, undifferentiated.
Declares the DissipatingBody class.
EvolModeType
The various evolution modes.
double semimajor() const
The current semimajor axis of the system.
const std::list< double > & semimajor_evolution_rate() const
The tabulated evolution of the semimajor axis so far.
Declares a class for a stopping condition that combines other stopping conditions.
std::valarray< Eigen::VectorXd > __above_lock_fractions_periapsis_deriv
The derivatives of the above lock fractions w.r.t. the periapses of the zones.
virtual void change_expansion_order(unsigned new_expansion_order)
Change the tidal potential expansion order for all dissipative zones.
A class combining the the outputs of multiple stopping conditions.
unsigned number_zones() const
The total number of zones in both system bodies.
double eccentricity() const
The current eccentricity of the system.
Describes a system of two bodies orbiting each other.