Planetary Orbital Evolution due to Tides
Orbital evolution of two objects experiencing tides
TidalPotentialExpansion.cpp
Go to the documentation of this file.
1 
9 
10 namespace Evolve {
11 
13  double radial_distance,
14  double azimuthal_angle,
15  double polar_angle,
16  double orbital_phase) const
17  {
18  double result = 0.0;
19  for(int m=-2; m<=2; ++m) {
20  std::complex<double> no_deriv,
21  inclination_deriv,
22  eccentricity_deriv,
23  error;
25  m,
26  mprime,
27  no_deriv,
28  inclination_deriv,
29  eccentricity_deriv);
30  assert(std::isfinite(no_deriv.real()));
31  result += (
32  no_deriv
33  *
34  std::pow(radial_distance, 2)
35  *
36  boost::math::spherical_harmonic(
37  2,
38  m,
39  polar_angle,
40  azimuthal_angle
41  )
42  *
43  std::complex<double>(std::cos(mprime * orbital_phase),
44  -std::sin(mprime * orbital_phase))
45  ).real();
46  assert(std::isfinite(result));
47  }
48  return result;
49  }
50 
52  double radial_distance,
53  double azimuthal_angle,
54  double polar_angle,
55  double time,
56  int expansion_order
57  )
58  {
59  std::cerr << "Approximating tidal potential to "
60  << expansion_order
61  << "order"
62  << std::endl;
63  assert(radial_distance >= 0);
64  assert(azimuthal_angle >= 0);
65  assert(azimuthal_angle < 2 * M_PI);
66  assert(polar_angle >= 0);
67  assert(polar_angle <= M_PI);
68 
70 
71  double orbital_phase = (
72  2.0 * M_PI * time
73  /
78  );
79 
80  double potential_norm = -(
82  *
84  /
85  std::pow(__semimajor, 3)
87 
88  double result = 0.0;
89  for(
90  int mprime = -expansion_order;
91  mprime <= expansion_order;
92  ++mprime
93  ) {
94  result += tidal_term(mprime,
95  radial_distance,
96  azimuthal_angle,
97  polar_angle,
98  orbital_phase);
99  }
100  return potential_norm * result;
101  }
102 
103 }//End Evolve namespace
Basic description of two bodies in an eccentric orbit.
double __semimajor
The semimajor axis of the orbit in solar radii.
double __eccentricity
The eccentricity of the orbit.
Declare an interface for evaluating the expansion of the tidal potential.
Orientations of zones of bodies in a binary system.
const double solar_radius
Radius of the sun [m].
void configure(double inclination, double arg_of_periapsis=0)
Set the inclination relative to the orbit.
TidalPotentialTerms __expansion_coef
The coefficients of the expansion of the tidal potential. ( )
const double solar_mass
Mass of the sun [kg].
double __secondary_mass
The mass of the perturber object in solar masses.
double tidal_term(int mprime, double radial_distance, double azimuthal_angle, double polar_angle, double orbital_phase) const
Return a single tidal term: .
double evaluate_spherical_coords(double radial_distance, double azimuthal_angle, double polar_angle, double time, int expansion_order)
double __primary_mass
The mass of the tidally perturbed object in solar masses.
double orbital_period() const
The orbital period of the system in days.
const double G
Gravitational constant in SI.