orbit.c ¶
Orbital evolution routines.
Functions
- 
void on_progress(progress_data_t *data)¶
- Update the progress bar. Used internally. 
- 
double modulus(double x, double y)¶
- The arithmetic modulus, x mod y. - Return
- x mod y
- Parameters
- x: A double
- y: A double
 
 
- 
double sgn(double x)¶
- Returns the sign of x. - Return
- sign(x)
- Parameters
- x: A double
 
 
- 
double TrueAnomaly(double E, double ecc)¶
- The true anomaly as a function of the eccentric anomaly E and the eccentricity ecc. - Return
- The true anomaly in radians
- Parameters
- E: The eccentric anomaly in radians
- ecc: The eccentricity
 
 
- 
double EccentricAnomalyFast(double dMeanA, double dEcc, double tol, int maxiter)¶
- Computes the eccentric anomaly from the mean anomaly and the eccentricity using a fast iterative scheme. Adapted from Rory Barnes, based on equations in Murray & Dermott. - Return
- The eccentric anomaly in radians
- Parameters
- dMeanA: The mean anomaly in radians
- dEcc: The eccentricity
- tol: The solver tolerance
- maxiter: The maximum number of iterations
 
 
- 
double EccentricAnomaly(double M, double e, double tol, int maxiter)¶
- Computes the eccentric anomaly from the mean anomaly and the eccentricity. This is a simpler version of the Kepler solver, borrowed from https://github.com/lkreidberg/batman/blob/master/c_src/_rsky.c - Return
- The eccentric anomaly in radians
- Parameters
- M: The mean anomaly in radians
- e: The eccentricity
- tol: The solver tolerance
- maxiter: The maximum number of iterations
 
 
- 
int Kepler(int np, BODY **body, SETTINGS settings)¶
- Compute the instantaneous x, y, and z positions of all planets with a simple multi-Keplerian solver. 
- 
void heartbeat(struct reb_simulation *r)¶
- Called at each step of the N-Body simulation. Currently does nothing! 
- 
int NBody(int np, BODY **body, SETTINGS settings, int halt_on_occultation, int occulted, int noccultors, int occultors[noccultors], int noccultations, double occultation_times[noccultations], int occultation_inds[noccultations], double occultation_durs[noccultations])¶
- Compute the instantaneous x, y, and z positions of all planets using the REBOUND N-Body code. - Return
- The error code
- Parameters
- np: The number of particles (bodies) in the system
- body: An array of BODY pointers to each of the bodies in the system
- settings: An instance of the SETTINGS class
- halt_on_occultation: Halt when a certain number of occultations of a given body have occurred?
- occulted: The index of the body being occulted (only used if- halt_on_occultation= 1)
- noccultors: Size of occultor array (only used if- halt_on_occultation= 1)
- occultors: The indices of the occultors of- occulted(only used if- halt_on_occultation= 1)
- noccultations: The number of occultations to look for (only used if- halt_on_occultation= 1)
- occultation_times: The times of each of the occultations, computed by this function (only used if- halt_on_occultation= 1)
- occultation_inds: The indices of the occultors corresponding to each of the occultations identified by this function (only used if- halt_on_occultation= 1)
- occultation_durs: The duration of each of the occultations identified by this function (only used if- halt_on_occultation= 1)