LM10Potential

class gala.potential.potential.LM10Potential(units=<UnitSystem (kpc, Myr, solMass, rad)>, disk={}, bulge={}, halo={})[source]

Bases: gala.potential.potential.ccompositepotential.CCompositePotential

The Galactic potential used by Law and Majewski (2010) to represent the Milky Way as a three-component sum of disk, bulge, and halo.

The disk potential is an axisymmetric MiyamotoNagaiPotential, the bulge potential is a spherical HernquistPotential, and the halo potential is a triaxial LogarithmicPotential.

Default parameters are fixed to those found in LM10 by fitting N-body simulations to the Sagittarius stream.

Parameters
unitsUnitSystem (optional)

Set of non-reducable units that specify (at minimum) the length, mass, time, and angle units.

diskdict (optional)

Parameters to be passed to the MiyamotoNagaiPotential.

bulgedict (optional)

Parameters to be passed to the HernquistPotential.

halodict (optional)

Parameters to be passed to the LogarithmicPotential.

Note: in subclassing, order of arguments must match order of potential
components added at bottom of init.

Attributes Summary

Wrapper

ndim

parameters

units

Methods Summary

__call__(q)

Call self as a function.

acceleration(q[, t])

Compute the acceleration due to the potential at the given position(s).

circular_velocity(q[, t])

Estimate the circular velocity at the given position assuming the potential is spherical.

clear()

copy()

density(q[, t])

Compute the density value at the given position(s).

energy(q[, t])

Compute the potential energy at the given position(s).

fromkeys(/, iterable[, value])

Create a new ordered dictionary with keys from iterable and values set to value.

get(key[, default])

Return the value for key if key is in the dictionary, else default.

gradient(q[, t])

Compute the gradient of the potential at the given position(s).

hessian(q[, t])

Compute the Hessian of the potential at the given position(s).

integrate_orbit(*args, **kwargs)

Warning

This is now deprecated. Convenient orbit integration should

items()

keys()

mass_enclosed(q, t)

Estimate the mass enclosed within the given position by assuming the potential is spherical.

move_to_end(/, key[, last])

Move an existing element to the end (or beginning if last is false).

plot_contours(grid[, filled, ax, labels, …])

Plot equipotentials contours.

plot_density_contours(grid[, filled, ax, …])

Plot density contours.

pop(k[,d])

value.

popitem(/[, last])

Remove and return a (key, value) pair from the dictionary.

replace_units(units)

Change the unit system of this potential.

save(f)

Save the potential to a text file.

setdefault(/, key[, default])

Insert key with a value of default if key is not in the dictionary.

to_galpy_potential([ro, vo])

Convert a Gala potential to a Galpy potential instance

to_latex()

Return a string LaTeX representation of this potential

to_sympy()

Return a representation of this potential class as a sympy expression

total_energy(x, v)

Compute the total energy (per unit mass) of a point in phase-space in this potential.

update([E, ]**F)

If E is present and has a .keys() method, then does: for k in E: D[k] = E[k] If E is present and lacks a .keys() method, then does: for k, v in E: D[k] = v In either case, this is followed by: for k in F: D[k] = F[k]

value(*args, **kwargs)

values()

Attributes Documentation

Wrapper = None
ndim = 3
parameters
units

Methods Documentation

__call__(q)

Call self as a function.

acceleration(q, t=0.0)

Compute the acceleration due to the potential at the given position(s).

Parameters
qPhaseSpacePosition, Quantity, array_like

Position to compute the acceleration at.

Returns
accQuantity

The acceleration. Will have the same shape as the input position array, q.

circular_velocity(q, t=0.0)

Estimate the circular velocity at the given position assuming the potential is spherical.

Parameters
qarray_like, numeric

Position(s) to estimate the circular velocity.

Returns
vcircQuantity

Circular velocity at the given position(s). If the input position has shape q.shape, the output energy will have shape q.shape[1:].

clear()None.  Remove all items from od.
copy()a shallow copy of od
density(q, t=0.0)

Compute the density value at the given position(s).

Parameters
qPhaseSpacePosition, Quantity, array_like

The position to compute the value of the potential. If the input position object has no units (i.e. is an ndarray), it is assumed to be in the same unit system as the potential.

Returns
densQuantity

The potential energy or value of the potential. If the input position has shape q.shape, the output energy will have shape q.shape[1:].

energy(q, t=0.0)

Compute the potential energy at the given position(s).

Parameters
qPhaseSpacePosition, Quantity, array_like

The position to compute the value of the potential. If the input position object has no units (i.e. is an ndarray), it is assumed to be in the same unit system as the potential.

Returns
EQuantity

The potential energy per unit mass or value of the potential.

fromkeys(/, iterable, value=None)

Create a new ordered dictionary with keys from iterable and values set to value.

get(key, default=None, /)

Return the value for key if key is in the dictionary, else default.

gradient(q, t=0.0)

Compute the gradient of the potential at the given position(s).

Parameters
qPhaseSpacePosition, Quantity, array_like

The position to compute the value of the potential. If the input position object has no units (i.e. is an ndarray), it is assumed to be in the same unit system as the potential.

Returns
gradQuantity

The gradient of the potential. Will have the same shape as the input position.

hessian(q, t=0.0)

Compute the Hessian of the potential at the given position(s).

Parameters
qPhaseSpacePosition, Quantity, array_like

The position to compute the value of the potential. If the input position object has no units (i.e. is an ndarray), it is assumed to be in the same unit system as the potential.

Returns
hessQuantity

The Hessian matrix of second derivatives of the potential. If the input position has shape q.shape, the output energy will have shape (q.shape[0],q.shape[0]) + q.shape[1:]. That is, an n_dim by n_dim array (matrix) for each position.

integrate_orbit(*args, **kwargs)

Warning

This is now deprecated. Convenient orbit integration should happen using the gala.potential.Hamiltonian class. With a static reference frame, you just need to pass your potential in to the Hamiltonian constructor.

Integrate an orbit in the current potential using the integrator class provided. Uses same time specification as Integrator.run() – see the documentation for gala.integrate for more information.

Parameters
w0PhaseSpacePosition, array_like

Initial conditions.

IntegratorIntegrator (optional)

Integrator class to use.

Integrator_kwargsdict (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.

**time_spec

Specification of how long to integrate. See documentation for parse_time_specification.

Returns
orbitOrbit
items()a set-like object providing a view on D’s items
keys()a set-like object providing a view on D’s keys
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.

move_to_end(/, key, last=True)

Move an existing element to the end (or beginning if last is false).

Raise KeyError if the element does not exist.

plot_contours(grid, filled=True, ax=None, labels=None, subplots_kw={}, **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
gridtuple

Coordinate grids or slice value for each dimension. Should be a tuple of 1D arrays or numbers.

filledbool (optional)

Use contourf() instead of contour(). Default is True.

axmatplotlib.Axes (optional)
labelsiterable (optional)

List of axis labels.

subplots_kwdict

kwargs passed to matplotlib’s subplots() function if an axes object is not specified.

kwargsdict

kwargs passed to either contourf() or plot().

Returns
figFigure
plot_density_contours(grid, filled=True, ax=None, labels=None, subplots_kw={}, **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
gridtuple

Coordinate grids or slice value for each dimension. Should be a tuple of 1D arrays or numbers.

filledbool (optional)

Use contourf() instead of contour(). Default is True.

axmatplotlib.Axes (optional)
labelsiterable (optional)

List of axis labels.

subplots_kwdict

kwargs passed to matplotlib’s subplots() function if an axes object is not specified.

kwargsdict

kwargs passed to either contourf() or plot().

Returns
figFigure
pop(k[, d])v, remove specified key and return the corresponding

value. If key is not found, d is returned if given, otherwise KeyError is raised.

popitem(/, last=True)

Remove and return a (key, value) pair from the dictionary.

Pairs are returned in LIFO order if last is true or FIFO order if false.

replace_units(units)

Change the unit system of this potential.

Parameters
unitsUnitSystem

Set of non-reducable units that specify (at minimum) the length, mass, time, and angle units.

save(f)

Save the potential to a text file. See save() for more information.

Parameters
fstr, file_like

A filename or file-like object to write the input potential object to.

setdefault(/, key, default=None)

Insert key with a value of default if key is not in the dictionary.

Return the value for key if key is in the dictionary, else default.

to_galpy_potential(ro=None, vo=None)

Convert a Gala potential to a Galpy potential instance

Parameters
roquantity-like (optional)
voquantity-like (optional)
classmethod to_latex()

Return a string LaTeX representation of this potential

Returns
latex_strstr

The latex expression as a Python string.

classmethod to_sympy()

Return a representation of this potential class as a sympy expression

Returns
exprsympy expression
varsdict

A dictionary of sympy symbols used in the expression.

total_energy(x, v)

Compute the total energy (per unit mass) of a point in phase-space in this potential. Assumes the last axis of the input position / velocity is the dimension axis, e.g., for 100 points in 3-space, the arrays should have shape (100, 3).

Parameters
xarray_like, numeric

Position.

varray_like, numeric

Velocity.

update([E, ]**F)None.  Update D from dict/iterable E and F.

If E is present and has a .keys() method, then does: for k in E: D[k] = E[k] If E is present and lacks a .keys() method, then does: for k, v in E: D[k] = v In either case, this is followed by: for k in F: D[k] = F[k]

value(*args, **kwargs)
values()an object providing a view on D’s values