BaseSphericalDensityProfile#

class profiles.density.BaseSphericalDensityProfile(**kwargs)[source]#

Abstract base class defining the shared behaviors of all spherically symmetric density profiles.

This class provides a unified interface for constructing and analyzing spherical density profiles common in astrophysics, such as those describing dark matter halos, stellar distributions, or gas profiles.

Several additional features are implemented to support both symbolic and numerical computations within a consistent, unit-aware framework.

Features#

Derived Profiles

The following symbolic derived profiles are automatically constructed and accessible via get_derived_profile():

  • derivative : Radial derivative of the density profile (from parent class).

  • enclosed_mass : The mass enclosed within radius \(r\).

    \[M(<r) = \int_0^r dr' 4\pi r'^2 \rho(r')\]
  • gravitational_field : Symbolic gravitational field at radius r.

    \[\nabla \Phi = \frac{-GM(<r)}{r^2} = \frac{-G}{r^2} \int_0^r dr' 4\pi r'^2 \rho(r')\]
  • gravitational_potential : Symbolic gravitational potential at radius r.

    \[\Phi = \frac{GM(<r)}{r} = \frac{G}{r} \int_0^r dr' 4\pi r'^2 \rho(r')\]
Numerical Operations

This class provides several numerical methods for common astrophysical computations:

Notes

  • This class is abstract; concrete subclasses must set __IS_ABSTRACT__ = False.

  • All subclasses define __function__ and __function_units__ to specify the symbolic density expression and its units.

  • Units are handled consistently using unyt.

  • Cosmological methods require astropy.cosmology and a valid cosmology.

  • Derived profiles follow the standard symbolic setup of BaseProfile.

Examples

Subclassing Example

from pisces.profiles.density import (
    BaseSphericalDensityProfile,
)
import sympy as sp
import unyt

class SimplePowerLawDensityProfile(
    BaseSphericalDensityProfile
):
    __IS_ABSTRACT__ = False
    __PARAMETERS__ = {
        "rho_0": 1.0 * unyt.Msun / unyt.kpc**3
    }

    @classmethod
    def __function__(cls, r, rho_0=1.0):
        return rho_0 / r**2

    @classmethod
    def __function_units__(
        cls,
        r,
        rho_0=unyt_quantity("1 Msun/kpc**3"),
    ):
        r, rho_0 = unyt.Unit(r), unyt.Unit(rho_0)
        return rho_0 / r**2

Numerical Usage

profile = SimplePowerLawDensityProfile(
    rho_0=0.1 * unyt.Msun / unyt.kpc**3
)

# Density at radius 1 kpc
rho = profile(1.0 * unyt.kpc)

# Enclosed mass at 5 kpc
M_enc = profile.compute_enclosed_mass(
    5 * unyt.kpc
)

# Surface density profile
Sigma = profile.compute_surface_density(
    2.0 * unyt.kpc
)

# Gravitational potential as symbolic derived profile
phi_profile = profile.get_derived_profile(
    "gravitational_potential"
)
phi_expr = phi_profile.__function__(
    sp.Symbol("r"), G=1.0, rho_0=1.0
)

Methods

__init__(**kwargs)

Initialize a profile instance with specific parameter values.

compute_circular_velocity(r[, units, G])

Compute the circular velocity at radius \(r\).

compute_cosmological_overdensity_profile(z, R)

Compute the spherical overdensity profile relative to the critical density at redshift \(z\).

compute_cosmological_overdensity_radius(z, ...)

Find the radius enclosing a target overdensity relative to the critical density at redshift \(z\).

compute_deflection_angle(R[, mode, units, ...])

Compute the gravitational lensing deflection angle at projected radius \(R\).

compute_einstein_radius(z_lens, z_source[, ...])

Compute the Einstein ring angular radius \(\\theta_E\) for a perfectly aligned source-lens system.

compute_enclosed_mass(r[, units])

Numerically compute the enclosed mass \(M(r)\) within radius \(r\).

compute_escape_velocity(r[, units, G])

Compute the escape velocity at radius \(r\).

compute_fractional_mass_radius(rmin, rmax[, ...])

Find the radius enclosing a given fraction of the total mass.

compute_gravitational_field(r[, units, G])

Numerically compute the gravitational field \(g(r)\) at radius \(r\).

compute_gravitational_potential(r[, units, G])

Numerically compute the gravitational potential \(\Phi(r)\).

compute_lensing_convergence(R, z_lens, z_source)

Compute the lensing convergence \(\kappa(R)\) at projected radius \(R\).

compute_surface_density(R[, units])

Numerically compute the projected surface density at radius R from the origin.

compute_total_mass([units])

Numerically compute the total mass of the profile by integrating to infinity.

from_dict(data)

Reconstruct a profile instance from a dictionary.

from_hdf5(h5obj[, name])

Reconstruct a profile from HDF5 attributes.

from_json(filepath)

Reconstruct a profile from a JSON file.

from_yaml(filepath)

Reconstruct a profile from a YAML file.

get_derived_profile(profile_name, **kwargs)

Access and instantiate a derived profile by name.

get_expression_latex([substitute])

Return the LaTeX representation of the profile's symbolic expression.

get_output_units(*argu)

Determine the output units of the operation given some set of input units.

get_parameters_latex()

Return a LaTeX table of the profile parameters.

lambdify_expression(expression)

Convert a symbolic expression into a callable function.

list_derived_profiles()

List all available derived profiles for this instance.

substitute_expression(expression)

Replace symbolic parameters with numerical values in an expression.

to_dict()

Serialize this profile to a minimal dictionary representation.

to_hdf5(h5obj[, name])

Store profile metadata into an HDF5 object as attributes.

to_json(filepath, **kwargs)

Serialize the profile to a JSON file.

to_yaml(filepath, **kwargs)

Serialize the profile to a YAML file.

Attributes

derived_profile_classes

Get the available derived profile classes for this instance.

parameter_symbols

Get the symbolic representations of the coordinate system parameters.

parameters

The parameters of this coordinate system.

variable_symbols

The symbols representing each of the coordinate axes in this coordinate system.

variables

The axes names present in this coordinate system.