BayesicFitting

Model Fitting and Evidence Calculation

View project on GitHub



class Kepplers2ndLaw( object )Source

Class for calculating Kepplers second law for planetary motion.

The projection of the orbit on the sky is not included in this class.

The algorithm was taken from Boule.

param abbr name limits comment
0 e eccentricity of orbit 0<e<1 0: circular
1 a semi major axis a>0
2 P period of the orbit P>0
3 T phase since periastron 0<T<2π

The parameters are initialized at [0.0, 1.0, 1.0, 0.0].

meanAnomaly( xdata, params )

Return the mean anomaly.

P = params[2] = period T = params[3] = periastron passage

M = 2 * pi * xdata / P - T

Parameters

  • xdata : array_like
         times in the orbit
  • params : array_like
         parameters: eccentr, semimajor, period, ppass

dMdx( xdata, params )
Return derivatives of M (mean anomaly) to xdata

Returns

  • dMdx : array_like
         derivatives of M to x (xdata)

dMdpar( xdata, params )
Return derivatives of M (mean anomaly) to relevant parameters.

Returns

  • dMdP : array_like
         derivatives of M to P (period)
  • dMdp : array_like
         derivatives of M to p (phase of periastron)

eccentricAnomaly( xdata, params, Estart=None )
Take the best one : Halleys method

It converges in a few iterations for e <= 0.999999999

eccentricAnomaly0( xdata, params )
Return the eccentric anomaly, i.e. the solution for E of

Standard method by Jean Meuss
  Astronomical Algorithms, 2nd ed.,
  Willmann-Bell, Inc, Virginia, 193-196, 397-399

e = params[0] = eccentricity M = mean anomaly

E = M + e * sin( E )

Parameters

  • xdata : array_like
         times in the orbit
  • params : array_like
         parameters: eccentr, semimajor, period, ppass

eccentricAnomaly1( xdata, params )
Newtons method.

Return the eccentric anomaly, i.e. the solution for E of

e = params[0] = eccentricity M = mean anomaly

E = M + e * sin( E )

Parameters

  • xdata : array_like
         times in the orbit
  • params : array_like
         parameters: eccentr, semimajor, period, ppass

eccentricAnomaly2( xdata, params, Estart=None )
Halleys method.

Return the eccentric anomaly, i.e. the solution for E of

e = params[0] = eccentricity M = mean anomaly

E = M + e * sin( E )

Parameters

  • xdata : array_like
         times in the orbit
  • params : array_like
         parameters: eccentr, semimajor, period, phase
  • Estart : array_like
         starting values for E

dEdM( xdata, params, cosE )
Return derivatives of E (eccentric anomaly) to mean anomaly

Returns

  • dEdM : array_like
         derivatives of E to M (mean anomaly)

dEdx( xdata, params, cosE )
Return derivatives of E (eccentric anomaly) to xdata

Returns

  • dEdx : array_like
         derivatives of E to x (xdata)

dEdpar( xdata, params, cosE, sinE )
Return derivatives of E (eccentric anomaly) to relevant parameters.

Returns

  • dEde : array_like
         derivatives of E to e (eccentricity)
  • dEdP : array_like
         derivatives of E to P (period)
  • dEdp : array_like
         derivatives of E to p (phase of periastron)

radiusAndTrueAnomaly( xdata, params )
Return the radius and the true anomaly.

e = params[0] = eccentricity a = params[1] = semimajor axis E = eccentric anomaly

r = a * ( 1 - e * cos( E ) )

v = 2 * arctan( sqrt( (1+e)/(1-e) ) * tan( E / 2 ) )

from Wikepedia => Trigoniometic Identities tan( E / 2 ) = sqrt( ( 1 - cos( E ) ) / ( 1 + cos( E ) ) )
              = sqrt( ( 1 - c ) * ( 1 + c ) / ( 1 + c )2 )
              = sqrt( s2 / ( 1 + c )2 )
              = s / ( 1 + c )
              = sin( E ) / ( 1 + cos( E ) ) Avoid cases where cos( E ) is too close to -1

Parameters

  • xdata : array_like
         times in the orbit
  • params : array_like
         parameters: eccentr, semimajor, inclin, ascpos, asclon, period, ppass

Returns

  • r : array_like
         radius
  • v : array_like
         true anomaly

drvdE( xdata, params, cosE, sinE )
Return derivatives of r (radius) and v (true anomaly) to eccentric anomaly

Parameters

  • xdata : array_like
         times in the orbit
  • params : array_like
         parameters: eccentr, semimajor, period, ppass
  • cosE : array_like
         cosine of E
  • sinE : array_like
         sine of E

Returns

  • drdE : array_like
         derivatives of r to E (eccentric anomaly)
  • dvdE : array_like
         derivatives of v to E (eccentric anomaly)

drvdx( xdata, params, cosE, sinE )
Return derivatives of r (radius) and v (true anomaly) to xdata

Returns

  • drdx : array_like
         derivatives of r to x (xdata)
  • dvdx : array_like
         derivatives of v to x (xdata)

drvdpar( xdata, params, E, cosE, sinE )
Return derivatives of r (radius) and v (true anomaly) to relevant parameters.

Returns

  • drde : array_like
         derivatives of r to e (eccentricity)
  • drda : array_like
         derivatives of r to a (semimajor axis)
  • drdP : array_like
         derivatives of r to P (period)
  • drdp : array_like
         derivatives of r to p (phase of periastron)
  • dvde : array_like
         derivatives of v to e (eccentricity)
  • dvdP : array_like
         derivatives of v to P (period)
  • dvdp : array_like
         derivatives of v to p (phase of periastron)

getMsini( stellarmass )
Return the mass of the exoplanet in Jupiter masses.

Parameters

  • stellarmass : float
         mass of the host star in solar masses.