Planetary Orbital Evolution due to Tides
Orbital evolution of two objects experiencing tides
EccentricOrbit.cpp
Go to the documentation of this file.
1 
8 #include "EccentricOrbit.h"
9 
10 double eccentric_anomaly_equation(double eccentric_anomaly,
11  void *orbital_phase_eccentricity)
12 {
13  double orbital_phase = reinterpret_cast<double*>(
14  orbital_phase_eccentricity
15  )[0];
16  double eccentricity = reinterpret_cast<double*>(
17  orbital_phase_eccentricity
18  )[1];
19  return (eccentric_anomaly
20  -
21  eccentricity * std::sin(eccentric_anomaly)
22  -
23  orbital_phase);
24 }
25 
26 namespace Evolve {
28  {
29  return (
31  /
32  (__primary_mass + __secondary_mass)
33  );
34  }
35 
36  double EccentricOrbit::eccentric_anomaly(double orbital_phase) const
37  {
38  const gsl_root_fsolver_type *FsolverType;
39  gsl_root_fsolver *solver;
40  double root = 0;
41  double lower_bound = 0.0, upper_bound = 2.0 * M_PI;
42  gsl_function Function;
43  orbital_phase -= 2.0 * M_PI * std::floor(orbital_phase / (2.0 * M_PI));
44  double params[] = {orbital_phase, __eccentricity};
45 
46  Function.function = &eccentric_anomaly_equation;
47  Function.params = params;
48 
49  FsolverType = gsl_root_fsolver_brent;
50  solver = gsl_root_fsolver_alloc (FsolverType);
51  gsl_root_fsolver_set(solver, &Function, lower_bound, upper_bound);
52 
53 #ifdef VERBOSE_DEBUG
54  std::cerr << std::setw(25) << "Iteartion"
55  << std::setw(25) << "LowerBound"
56  << std::setw(25) << "root"
57  << std::setw(25) << "UpperBound"
58  << std::endl;
59 #endif
60 
61  int status = GSL_CONTINUE;
62  for(unsigned iter = 0; status == GSL_CONTINUE; ++iter) {
63  status = gsl_root_fsolver_iterate(solver);
64  root = gsl_root_fsolver_root(solver);
65  lower_bound = gsl_root_fsolver_x_lower(solver);
66  upper_bound = gsl_root_fsolver_x_upper(solver);
67  status = gsl_root_test_interval(lower_bound,
68  upper_bound,
69  1e-12,
70  1e-12);
71 #ifdef VERBOSE_DEBUG
72  std::cerr << std::setw(25) << iter
73  << std::setw(25) << lower_bound
74  << std::setw(25) << root
75  << std::setw(25) << upper_bound
76  << std::endl;
77 #endif
78  }
79 
80  gsl_root_fsolver_free(solver);
81 
82  assert(status == GSL_SUCCESS);
83 
84  return root;
85  }
86 
87  Eigen::Matrix<long double, 3, 1> EccentricOrbit::secondary_position(
88  double orbital_phase
89  ) const
90  {
91  long double current_eccentric_anomaly = eccentric_anomaly(orbital_phase);
92  long double semimajor = __semimajor,
93  eccentricity = __eccentricity,
94  one = 1.0;
95  typedef long double ldbl;
96  return Eigen::Matrix<ldbl, 3, 1>(
97  semimajor * (std::cos(current_eccentric_anomaly)
98  -
99  eccentricity),
100  semimajor * (sqrt(one - std::pow(eccentricity, 2))
101  *
102  std::sin(current_eccentric_anomaly)),
103  0.0
104  );
105  }
106 
108  {
109  return (
111  *
113  *
114  (2.0 * M_PI / (orbital_period() * Core::AstroConst::day))
115  *
116  std::sqrt(1.0 - std::pow(__eccentricity, 2))
117  );
118  }
119 
121  {
122  return (
124  *
126  *
128  /
130  );
131  }
132 
134  {
135  return std::sqrt(
136  4.0 * M_PI * M_PI
137  *
139  /
140  (
142  *
144  *
146  )
148  }
149 }
double reduced_mass() const
The reduced mass of the system in solar masses.
double eccentricity() const
The semimajor axis of the system.
double semimajor() const
The semimajor axis of the system.
Orientations of zones of bodies in a binary system.
double orbital_energy() const
The magnitude of the orbital energy of the system in .
double __eccentricity
The eccentricity of the orbit.
const double solar_radius
Radius of the sun [m].
double eccentric_anomaly(double orbital_phase) const
double __primary_mass
The mass of the tidally perturbed object in solar masses.
const double day
day [s]
double orbital_angmom() const
The magnitude of the orbital angular momentum of the system in .
const double solar_mass
Mass of the sun [kg].
double __secondary_mass
The mass of the perturber object in solar masses.
Eigen::Matrix< long double, 3, 1 > secondary_position(double orbital_phase) const
Secondary position vector in a coordinate system centered on the primary, with and ...
double orbital_period() const
The orbital period of the system in days.
Declare an interface for working with eccentric orbits.
const double G
Gravitational constant in SI.
double __semimajor
The semimajor axis of the orbit in solar radii.