Orbit

class pygro.orbit.Orbit[source]View on GitHub

Base class for integrating physical orbits in spherically symmetric spacetimes. It is a wrapper around the Geodesic object for the time-like case with additional methods to assign initial conditions based on the standard Keplerian parameters in spherically symmetric space-times.

__init__(geo_engine: GeodesicEngine | None = None, verbose: bool = False, orbital_parameters: dict | None = None, radial_coordinate: Basic | None = None, latitude: Basic | None = None, longitude: Basic | None = None)[source]View on GitHub

The Orbit constructor accepts the following arguments:

Parameters:
  • geo_engine (Optional[GeodesicEngine]) – The GeodesicEngine object to link to the Orbit. If not provided, the Geodesic will be linked to the last initialized GeodesicEngine.

  • verbose (bool) – Specifies whether to log information on the geodesic status to the standard output.

  • orbital_parameters (Optional[dict]) – An optional dictionary to be passed to the set_orbital_parameters() method. See the documentation of this method for the exact details of what to include in the dictionary.

  • radial_coordinate (Optional[sp.Basic]) – When initializing the Orbit, it will assume by the default that the radial coordinate in the metric will be the first component of the coordinate array (i.e. the content of Metric.x[1]). This can be overridden by specifically passing the symbol corresponding to the radial coordinate.

  • latitude (Optional[sp.Basic]) – When initializing the Orbit, it will assume by the default that the latitude angular coordinate (\(\theta\) in the usual Schwarzschild coordinates, see Define your own space-time) in the metric will be the second component of the coordinate array (i.e. the content of Metric.x[2]). This can be overridden by specifically passing the symbol corresponding to the latitude.

  • longitude (Optional[sp.Basic]) – When initializing the Orbit, it will assume by the default that the longitude angular coordinate (\(\phi\) in the usual Schwarzschild coordinates, see Define your own space-time) in the metric will be the third component of the coordinate array (i.e. the content of Metric.x[3]). This can be overridden by specifically passing the symbol corresponding to the longitude.

After initialization, the Orbit instance has the following attributes:

Variables:
  • params (dict) – The dictionary containing the orbital parameters

  • V_eff (Callable) – A callable which returns the effective potential of the orbit as a function of the orbital energy, angular momentum and of the radial coordinate \(V_\textrm{eff} = V_\textrm{eff}(E, L, r)\).

  • E (float) – The value of the orbital energy of the orbit that corresponds to the given orbital parameters.

  • L (float) – The value of the orbital angular momentum of the orbit that corresponds to the given orbital parameters.

integrate(tauf: float | int, initial_step: float | int, verbose: bool = False, direction: Literal['bw', 'fw'] = 'fw', initial_conditions_tol: float = 1e-15, **integration_kwargs)[source]View on GitHub

A wrapper around the integrate() method.

It accepts the same arguments plus:

Parameters:

initial_conditions_tol (float) – sets the numerical tolerance for the root finding algorithm of the initial conditions.

This method sets the initial conditions of the integration starting from the orbital parameters. In case the procedure is not successfull it falls back to the classical Keplerian initial conditions (alerting the user that this has been the case).

set_orbital_parameters(t_P: float | None = None, a: float | None = None, e: float | None = None, i: float | None = None, omega: float | None = None, Omega: float | None = None, t_A: float | None = None)[source]View on GitHub

Sets the orbital parameters that uniquely identify the orbit on the equatorial plane.

For a detailed explaination on how a choice of orbital parameters is translated into intial conditions for the geodesic, see the Integrating physical orbits tutorial and the example notebooks Circular orbits and Orbital precession.

The parametrization that we choose is the usual one in celestial mechanics:

  • \(t_p\) (t_P): coordinate time of pericenter passage. Alternatively, one can specify the coordinate time of apocenter passge \(t_a\) (t_A). Assigning one (and only one) of these parameters is mandatory and will raise an error if not assigned. Depending on the choice between t_P and t_A the integration will be started at pericenter or apocenter.

  • \(a\) (a): The semi-major axis of the orbit in the same units in which the radial coordiante in the metrix is expressed. This parameter is mandatory and will raise an error if not assigned.

  • \(e\) (e): The orbital eccentricity, which fixes the radial turning points \(r_p=a(1-e)\) and \(r_a=a(1+e)\) (see Integrating physical orbits). This parameter is mandatory and will raise an error if not assigned.

  • \((i,\,\omega,\,\Omega)\) (i, omega, Omega): the three angular orbital parameters (in radians) which will rotate the initial conditions out of the equatorial plane. See Integrating physical orbits for a proper definition of these parameters. These parameters are optional. If not specified they will be fixed to zero, leaving the orbit on the equatorial plane (\(\theta = \pi/2\) or the corresponding latitude coordinate).