KuzminPotential#
- class gala.potential.potential.KuzminPotential(m, a, units=None, origin=None, R=None)[source]#
Bases:
CPotentialBaseKuzmin potential for a flattened mass distribution.
- Parameters:
- m
Quantity,numeric[mass] Mass.
- a
Quantity,numeric[length] Flattening parameter.
- units
UnitSystem(optional) Set of non-reducable units that specify (at minimum) the length, mass, time, and angle units.
- origin
Quantity(optional) The origin of the potential, the default being 0.
- R
Rotation, array_like (optional) A Scipy
Rotationobject or an array representing a rotation matrix that specifies a rotation of the potential. This is applied after the origin shift. Default is the identity matrix.
- m
Attributes Summary
Methods Summary
__call__(q)Call self as a function.
acceleration([q, t])Compute the gravitational acceleration at the given position(s).
as_interop(package, **kwargs)Interoperability with other Galactic dynamics packages
circular_velocity([q, t])Estimate the circular velocity at given position(s) assuming spherical symmetry.
density([q, t])Compute the mass density at the given position(s).
energy([q, t])Compute the gravitational potential energy at the given position(s).
gradient([q, t])Compute the gradient of the gravitational potential.
hessian([q, t])Compute the Hessian matrix of the gravitational potential.
integrate_orbit(*args, **kwargs)Integrate an orbit in the current potential using the integrator class provided.
mass_enclosed(q, t)Estimate the mass enclosed within the given position by assuming the potential is spherical.
plot_contours(grid[, t, filled, ax, labels, ...])Plot equipotentials contours.
plot_density_contours(grid[, t, filled, ax, ...])Plot density contours.
plot_rotation_curve(R_grid[, t, ax, labels])Plot the rotation curve or circular velocity curve for this potential on the input grid of cylindrical radii.
replace_units(units)Change the unit system of this potential.
replicate(**kwargs)Return a copy of the potential instance with some parameter values changed.
save(f)Save the potential to a text file.
to_latex()Return a string LaTeX representation of this potential
to_sympy(v, p)Return a representation of this potential class as a sympy expression
Attributes Documentation
- a = <PotentialParameter: a [length]>#
- m = <PotentialParameter: m [mass]>#
- ndim = 3#
- units#
Methods Documentation
- __call__(q)#
Call self as a function.
- acceleration(q=None, t=0.0, **coord_kwargs)#
Compute the gravitational acceleration at the given position(s).
The gravitational acceleration is the negative gradient of the potential: \(\vec{a} = -\nabla \phi\). This is the acceleration experienced by a test particle in the gravitational field.
- Parameters:
- q
PhaseSpacePosition,Quantity, array_like, optional Position(s) at which to compute the gravitational acceleration. If the input has no units (i.e., is an
ndarray), it is assumed to be in the same unit system as the potential. If using symmetry coordinates, pass those as keyword arguments instead and leave q as None.- t
numeric,Quantity, optional Time at which to evaluate the acceleration. Default is 0.
- **coord_kwargs
For potentials with spherical or cylindrical symmetry, you can optionally provide coordinates in the natural coordinate system. For spherical potentials, use
r=.... For cylindrical potentials, useR=...and optionallyz=...(defaults to 0).
- q
- Returns:
- acc
Quantity The gravitational acceleration vector(s). Has the same shape as the input position array
q. Units are acceleration (e.g., m/s² in SI units).
- acc
See also
gradientCompute the potential gradient (negative acceleration).
Notes
This method is equivalent to
-self.gradient(q, t)and is provided for convenience in orbital integration and dynamics calculations.
- as_interop(package, **kwargs)#
Interoperability with other Galactic dynamics packages
- Parameters:
- package
str The package to export the potential to. Currently supported packages are
"galpy"and"agama".- kwargs
Any additional keyword arguments are passed to the interop function.
- package
- circular_velocity(q=None, t=0.0, **coord_kwargs)#
Estimate the circular velocity at given position(s) assuming spherical symmetry.
The circular velocity is the speed required for a circular orbit at the given radius in a spherically symmetric potential. It is computed using \(v_{\rm circ}(r) = \sqrt{r |dΦ/dr|}\), where the radial derivative is evaluated from the potential gradient.
- Parameters:
- q
PhaseSpacePosition,Quantity, array_like, optional Position(s) at which to estimate the circular velocity. The calculation uses the spherical radius from the origin. If the input has no units, it is assumed to be in the same unit system as the potential. If using symmetry coordinates, pass those as keyword arguments instead and leave q as None.
- t
numeric,Quantity, optional Time at which to evaluate the circular velocity. Default is 0.
- **coord_kwargs
For potentials with spherical or cylindrical symmetry, you can optionally provide coordinates in the natural coordinate system. For spherical potentials, use
r=.... For cylindrical potentials, useR=...and optionallyz=...(defaults to 0).
- q
- Returns:
- vcirc
Quantity Circular velocity at the spherical radius corresponding to each position. For input shape
(n_dim, n_positions), returns shape(n_positions,). Units are velocity (e.g., m/s in SI units).
- vcirc
Notes
This method assumes the potential is approximately spherically symmetric. The circular velocity is computed using the relation:
\[v_{\rm circ}(r) = \sqrt{r \left| \frac{d\Phi}{dr} \right|}\]where the radial derivative is computed from the Cartesian gradient via \(dΦ/dr = \vec{\nabla}Φ \cdot \hat{r}\).
For exactly spherical potentials, this gives the speed of circular orbits. For non-spherical potentials, this provides an approximation useful for initial orbit estimates.
- density(q=None, t=0.0, **coord_kwargs)#
Compute the mass density at the given position(s).
For potentials that have an associated mass distribution, this method computes the mass density rho(q, t) at the specified positions and time. The density is related to the potential via Poisson’s equation: \(\nabla^2 \phi = 4\pi G \rho\).
- Parameters:
- q
PhaseSpacePosition,Quantity, array_like, optional Position(s) at which to evaluate the mass density. If the input has no units (i.e., is an
ndarray), it is assumed to be in the same unit system as the potential. Shape should be(n_dim,)for a single position or(n_dim, n_positions)for multiple positions. If using symmetry coordinates, pass those as keyword arguments instead and leave q as None.- t
numeric,Quantity, optional Time at which to evaluate the mass density. Default is 0.
- **coord_kwargs
For potentials with spherical or cylindrical symmetry, you can optionally provide coordinates in the natural coordinate system. For spherical potentials, use
r=.... For cylindrical potentials, useR=...and optionallyz=...(defaults to 0).
- q
- Returns:
- dens
Quantity The mass density at the specified position(s). For input shape
(n_dim, n_positions), returns shape(n_positions,). Units are mass density (e.g., kg/m³ in SI units).
- dens
- Raises:
NotImplementedErrorIf the potential does not have an implemented density function.
Notes
Not all potential models have an implemented density function. For potentials without a density implementation, this method will raise a
NotImplementedError.The density is computed using the relationship with the potential’s Laplacian (when available) or from the underlying mass model.
- energy(q=None, t=0.0, **coord_kwargs)#
Compute the gravitational potential energy at the given position(s).
The potential energy per unit mass is evaluated at the specified position(s) and time.
- Parameters:
- q
PhaseSpacePosition,Quantity, array_like, optional Position(s) at which to evaluate the potential. If the input has no units (i.e., is an
ndarray), it is assumed to be in the same unit system as the potential. Shape should be(n_dim,)for a single position or(n_dim, n_positions)for multiple positions. If using symmetry coordinates, pass those as keyword arguments instead and leave q as None.- t
numeric,Quantity, optional Time at which to evaluate the potential. Default is 0.
- **coord_kwargs
For potentials with spherical or cylindrical symmetry, you can optionally provide coordinates in the natural coordinate system. For spherical potentials, use
r=.... For cylindrical potentials, useR=...and optionallyz=...(defaults to 0).
- q
- Returns:
- E
Quantity The gravitational potential energy per unit mass. For input shape
(n_dim, n_positions), returns shape(n_positions,). Units are specific energy (e.g., m²/s² in SI units).
- E
Notes
The potential energy is related to the gravitational acceleration by \(\vec{a} = -\nabla \phi\), where φ is the potential energy per unit mass.
Examples
Using Cartesian coordinates (works for all potentials):
>>> import astropy.units as u >>> import numpy as np >>> pot = SomePotential(...) >>> xyz = np.array([[1., 0., 0.]]).T * u.kpc >>> pot.energy(xyz)
For spherical potentials, you can use spherical radius:
>>> pot = HernquistPotential(m=1e10*u.Msun, c=5*u.kpc) >>> r = np.linspace(0.1, 10, 100) * u.kpc >>> pot.energy(r=r)
For cylindrical potentials, you can use R and z:
>>> pot = MiyamotoNagaiPotential(m=1e11*u.Msun, a=3*u.kpc, b=0.3*u.kpc) >>> R = np.linspace(1, 15, 100) * u.kpc >>> pot.energy(R=R, z=0*u.kpc) >>> pot.energy(R=R) # z defaults to 0
- gradient(q=None, t=0.0, **coord_kwargs)#
Compute the gradient of the gravitational potential.
- Parameters:
- q
PhaseSpacePosition,Quantity, array_like, optional Position(s) at which to evaluate the potential gradient. If the input has no units (i.e., is an
ndarray), it is assumed to be in the same unit system as the potential. Shape should be(n_dim,)for a single position or(n_dim, n_positions)for multiple positions. If using symmetry coordinates, pass those as keyword arguments instead and leave q as None.- t
numeric,Quantity, optional Time at which to evaluate the potential gradient. Default is 0.
- **coord_kwargs
For potentials with spherical or cylindrical symmetry, you can optionally provide coordinates in the natural coordinate system. For spherical potentials, use
r=.... For cylindrical potentials, useR=...and optionallyz=...(defaults to 0).
- q
- Returns:
- grad
Quantity The gradient of the gravitational potential. Has the same shape as the input position array
q. Units are acceleration (e.g., m/s² in SI units). To get gravitational acceleration, useacceleration()or compute-gradient().
- grad
See also
accelerationCompute gravitational acceleration (negative gradient).
Notes
The relationship between potential φ, gradient, and acceleration is:
\[\vec{a} = -\nabla \phi = -\frac{\partial \phi}{\partial \vec{q}}\]The gradient is always returned in Cartesian coordinates, even when using symmetry coordinates as input.
- hessian(q=None, t=0.0, **coord_kwargs)#
Compute the Hessian matrix of the gravitational potential.
The Hessian matrix contains the second partial derivatives of the potential: \(H_{ij} = \frac{\partial^2 \phi}{\partial q_i \partial q_j}\). This is useful for stability analysis, computing tidal tensors, and orbital frequency analysis.
- Parameters:
- q
PhaseSpacePosition,Quantity, array_like, optional Position(s) at which to evaluate the Hessian matrix. If the input has no units (i.e., is an
ndarray), it is assumed to be in the same unit system as the potential. Shape should be(n_dim,)for a single position or(n_dim, n_positions)for multiple positions. If using symmetry coordinates, pass those as keyword arguments instead and leave q as None.- t
numeric,Quantity, optional Time at which to evaluate the Hessian matrix. Default is 0.
- **coord_kwargs
For potentials with spherical or cylindrical symmetry, you can optionally provide coordinates in the natural coordinate system. For spherical potentials, use
r=.... For cylindrical potentials, useR=...and optionallyz=...(defaults to 0).
- q
- Returns:
- hess
Quantity The Hessian matrix of second derivatives. For input shape
(n_dim, n_positions), returns shape(n_dim, n_dim, n_positions). Eachn_dim x n_dimslice corresponds to the Hessian matrix at one position. Units are acceleration per length (e.g., s⁻² in SI units).
- hess
- Raises:
NotImplementedErrorIf the potential does not have an implemented Hessian function, or if the potential is rotated (
Ris not the identity).
Notes
Computing Hessian matrices for rotated potentials (when
Ris not the identity matrix) is currently not supported and will raise aNotImplementedError.Not all potential models have an implemented Hessian function. For potentials without a Hessian implementation, this method will raise a
NotImplementedError.The Hessian matrix is symmetric for time-independent potentials.
- integrate_orbit(*args, **kwargs)#
Integrate an orbit in the current potential using the integrator class provided. Uses same time specification as
Integrator()– see the documentation forgala.integratefor more information.- Parameters:
- w0
PhaseSpacePosition, array_like Initial conditions.
- Integrator
Integrator(optional) Integrator class to use.
- Integrator_kwargs
dict(optional) Any extra keyword argumets to pass to the integrator class when initializing. Only works in non-Cython mode.
- cython_if_possiblebool (optional)
If there is a Cython version of the integrator implemented, and the potential object has a C instance, using Cython will be much faster.
- save_allbool (optional)
Controls whether to store the phase-space position at all intermediate timesteps. Set to False to store only the final values (i.e. the phase-space position(s) at the final timestep). Default is True.
- **time_spec
Specification of how long to integrate. See documentation for
parse_time_specification.
- w0
- Returns:
- orbit
Orbit
- orbit
- mass_enclosed(q, t)#
Estimate the mass enclosed within the given position by assuming the potential is spherical. This is not so good!
- Parameters:
- qarray_like,
numeric Position to compute the mass enclosed.
- qarray_like,
- plot_contours(grid, t=0.0, filled=True, ax=None, labels=None, subplots_kw=None, **kwargs)#
Plot equipotentials contours. Computes the potential energy on a grid (specified by the array
grid).Warning
Right now the grid input must be arrays and must already be in the unit system of the potential. Quantity support is coming…
- Parameters:
- grid
tuple Coordinate grids or slice value for each dimension. Should be a tuple of 1D arrays or numbers.
- tquantity-like (optional)
The time to evaluate at.
- filledbool (optional)
Use
contourf()instead ofcontour(). Default isTrue.- ax
matplotlib.Axes(optional) - labelsiterable (optional)
List of axis labels.
- subplots_kw
dict kwargs passed to matplotlib’s subplots() function if an axes object is not specified.
- kwargs
dict kwargs passed to either contourf() or plot().
- grid
- Returns:
- fig
Figure
- fig
- plot_density_contours(grid, t=0.0, filled=True, ax=None, labels=None, subplots_kw=None, **kwargs)#
Plot density contours. Computes the density on a grid (specified by the array
grid).Warning
For now, the grid input must be arrays and must already be in the unit system of the potential. Quantity support is coming…
- Parameters:
- grid
tuple Coordinate grids or slice value for each dimension. Should be a tuple of 1D arrays or numbers.
- tquantity-like (optional)
The time to evaluate at.
- filledbool (optional)
Use
contourf()instead ofcontour(). Default isTrue.- ax
matplotlib.Axes(optional) - labelsiterable (optional)
List of axis labels.
- subplots_kw
dict kwargs passed to matplotlib’s subplots() function if an axes object is not specified.
- kwargs
dict kwargs passed to either contourf() or plot().
- grid
- Returns:
- fig
Figure
- fig
- plot_rotation_curve(R_grid, t=0.0, ax=None, labels=None, **plot_kwargs)#
Plot the rotation curve or circular velocity curve for this potential on the input grid of cylindrical radii.
- Parameters:
- R_gridarray_like
A grid of radius values to compute the rotation curve at. This should be a one-dimensional grid.
- tquantity-like (optional)
The time to evaluate at.
- ax
matplotlib.Axes(optional) - labelsiterable (optional)
List of axis labels. Set to False to disable adding labels.
- plot_kwargs
dict kwargs passed to plot().
- Returns:
- fig
Figure - ax
Axes
- fig
- replace_units(units)#
Change the unit system of this potential.
- Parameters:
- units
UnitSystem Set of non-reducable units that specify (at minimum) the
- length, mass, time, and angle units.
- units
- replicate(**kwargs)#
Return a copy of the potential instance with some parameter values changed. This always produces copies of any parameter arrays.
- Parameters:
- **kwargs
All other keyword arguments are used to overwrite parameter values when making the copy.
- Returns:
- replicant
PotentialBasesubclass instance The replicated potential.
- replicant
- save(f)#
Save the potential to a text file. See
save()for more information.- Parameters:
- f
str,file_like A filename or file-like object to write the input potential object to.
- f