Unitless quantities methods

class pNbody._unitless_quantities_mixin._NbodyUnitlessQuantitiesMixin

Bases: object

Mixin providing methods to compute intrinsic physical properties from raw simulation data.

This class calculates quantities derived directly from the state vectors (positions, velocities, masses, etc.) stored in the simulation snapshots. These properties are ‘dimensionless’ or ‘unit-less,’ meaning they are expressed in the simulation’s internal coordinate system.

To obtain values in physical units (e.g., kpc, Msun, km/s), the outputs of these methods must be scaled by the appropriate transformation factors defined by the simulation’s normalization.

Accel(x, eps)

Compute the gravitational acceleration at specific target positions.

Calculates the acceleration vector $mathbf{a}(mathbf{x}) = sum

rac{m_i (mathbf{x}_i - mathbf{x})}{(|mathbf{x} - mathbf{x}_i|^2 + epsilon^2)^{3/2}}$

by direct summation of all particle contributions.

xndarray or list

Target position(s) where the acceleration is evaluated. Can be a single vector [x, y, z] or an array of shape (N, 3).

epsfloat

The softening length used to smooth the gravitational force at short ranges.

ndarray

The acceleration vector(s) at the input position(s). Returns shape (3,) for a single input or (N, 3) for multiple targets.

The acceleration is calculated in internal code units where the gravitational constant $G$ is assumed to be $1$. The calculation is parallelized via MPI across the particle distribution.

Ekin()

Calculate the total kinetic energy of the system.

Computes the sum of $1/2 m v^2$ for all particles across all MPI processes.

Returns

Total kinetic energy ($E_{kin}$) in internal units.

Return type

float

Epot(eps)

Calculate the total gravitational potential energy.

Computes the pairwise potential energy sum using a Plummer-like softening kernel to avoid numerical divergences at short distances.

Parameters

eps (float) – The softening length (gravitational smoothing scale).

Returns

Total potential energy ($W$ or $Phi$) in internal units.

Return type

float

Warning

This method is currently serial and does not support MPI parallelization. It should only be called on a single process.

L()

Calculate the angular momentum vector for each particle.

Computes $mathbf{L} = m (mathbf{r} imes mathbf{v})$ for every particle in the snapshot.

Returns

A (N, 3) float array containing the angular momentum components $[L_x, L_y, L_z]$ for each particle.

Return type

ndarray

Ltot()

Calculate the total angular momentum of the entire system.

Performs a global MPI reduction to sum the angular momentum vectors of all particles across all processes: $mathbf{L}_{tot} = sum mathbf{L}_i$.

Returns

A (3,) float array representing the total angular momentum vector $[L_{x,tot}, L_{y,tot}, L_{z,tot}]$.

Return type

ndarray

Pot(x, eps)

Compute the gravitational potential at one or more target positions.

Calculates the potential $Phi(mathbf{x}) = -sum

rac{m_i}{sqrt{|mathbf{x} - mathbf{x}_i|^2 + epsilon^2}}$

by summing contributions from all particles in the system.

xndarray or list

The target position(s) where the potential is evaluated. Can be a single vector [x, y, z] or an array of shape (N, 3).

epsfloat

The softening length used to smooth the gravitational force at short ranges.

float or ndarray

The potential value at the input position(s).

The potential is calculated in internal code units where the gravitational constant $G$ is assumed to be $1$. For physical units, multiply the result by the appropriate value of $G$.

property R

Alias for rxy().

TreeAccel(pos, eps, Tree=None)

Compute the gravitational acceleration using an Octree-based approximation.

Accelerates the calculation using a Barnes-Hut tree structure to approximate distant particle groups as monopoles (or higher orders).

Parameters
  • pos (ndarray) – A (N, 3) array of positions where the acceleration should be evaluated.

  • eps (float) – The softening length (gravitational smoothing scale).

  • Tree (Tree object, optional) – A pre-computed gravitational tree. If None, a new tree is constructed from the current particle distribution. Default is None.

Returns

An array of acceleration vectors [ax, ay, az] for each target position.

Return type

ndarray

Notes

The acceleration is calculated in internal code units where the gravitational constant $G$ is assumed to be $1$.

Warning

This method is currently serial and does not support MPI parallelization.

TreePot(pos, eps, Tree=None)

Compute the gravitational potential using an Octree-based approximation.

Uses a Barnes-Hut tree algorithm to estimate the potential at target positions, significantly reducing computational cost from $O(N^2)$ to $O(N log N)$.

Parameters
  • pos (ndarray) – A (N, 3) array of positions where the potential should be evaluated.

  • eps (float) – The softening length (gravitational smoothing scale).

  • Tree (Tree object, optional) – A pre-computed gravitational tree. If None, a new tree is constructed from the current particle distribution. Default is None.

Returns

An array of potential values at the requested positions.

Return type

ndarray

Notes

The potential is calculated in internal code units where the gravitational constant $G$ is assumed to be $1$.

Warning

This method is currently serial and does not support MPI parallelization.

VR()

Calculate the radial velocity in the cylindrical coordinate system.

Compute the component of the velocity vector along the projected radial vector in the xy-plane. This represents the rate at which the cylindrical radius $R$ is changing.

Notes

This method is a functional alias for Vrxy().

Returns

A 1D array of cylindrical radial velocities ($v_R$) in dimensionless internal units. Positive values indicate outward motion from the z-axis.

Return type

ndarray

VT()

Calculate the tangential velocity in the cylindrical coordinate system.

Compute the velocity component perpendicular to the projected radial vector in the xy-plane. This represents the circular or azimuthal motion around the z-axis.

Notes

This method is a functional alias for Vtxy(). Note that the sign convention (direction of rotation) is preserved from the original method.

Returns

A 1D array of cylindrical tangential velocities ($v_{phi}$) in dimensionless internal units.

Return type

ndarray

Vphi()

Calculate the azimuthal velocity component (v_phi) in spherical coordinates.

Compute the velocity component in the direction of increasing longitude (azimuth). This represents the circular motion of particles around the z-axis, projected into the xy-plane.

Returns

A 1D array of azimuthal velocities ($v_{phi}$) in dimensionless internal units. Positive values indicate counter-clockwise rotation around the z-axis.

Return type

ndarray

Vr()

Calculate the radial velocity in the spherical coordinate system.

Compute the component of the velocity vector along the radial line from the origin to each particle. This represents the rate of change of the spherical distance $r$. This is an alias for vrxyz().

Notes

This method is a functional alias for vrxyz().

Returns

A 1D array of spherical radial velocities ($v_r$) in dimensionless internal units. Positive values indicate outward radial motion.

Return type

ndarray

Vrxy()

Calculate the radial velocity in the cylindrical coordinate system.

Compute the projection of the velocity vector onto the radial unit vector in the xy-plane (z=0). This represents the expansion or contraction rate of the cylindrical radius $R$.

Returns

A 1D array of cylindrical radial velocities ($v_R$) in dimensionless internal units. Positive values indicate motion away from the z-axis.

Return type

ndarray

Vt()

Calculate the tengential velocity in the spherical coordinate system.

Compute the velocity component perpendicular to the radial vector. In 3D spherical coordinates, this represents the combined motion from both the polar (theta) and azimuthal (phi) angular changes.

It is calculated using the Pythagorean identity: $v_{tangential} = sqrt{v_{total}^2 - v_{radial}^2}$

Returns

A 1D array of tangential velocity magnitudes ($v_t$) in dimensionless internal units.

Return type

ndarray

Vtheta()

Calculate the polar velocity component (v_theta) in spherical coordinates.

Compute the velocity component in the direction of increasing polar angle (theta), measured down from the positive z-axis (co-latitude). This represents the motion of particles along lines of longitude.

Returns

A 1D array of polar velocities ($v_{ heta}$) in dimensionless internal units. Positive values indicate motion moving away from the North pole (increasing theta).

Return type

ndarray

Vtxy()

Calculate the tangential (circular) velocity in the cylindrical coordinate system.

Compute the component of the velocity vector perpendicular to the radial vector in the cylindrical coordinate system ($z=0$ plane). This represents the instantaneous rotational speed around the z-axis.

Returns

A 1D array of tangential velocities ($v_{phi}$) in dimensionless internal units. Positive values indicate counter-clockwise rotation.

Return type

ndarray

Vx()

Return the x-components of the velocity.

Extract the first column of the velocity array.

Notes

This method is a functional alias for vx().

Returns

A 1D array of x-velocities in dimensionless internal units.

Return type

ndarray

Vy()

Return the y-components of the velocity.

Extract the second column of the velocity array.

Notes

This method is a functional alias for vy().

Returns

A 1D array of y-velocities in dimensionless internal units.

Return type

ndarray

Vz()

Return the z-components of the velocity.

Extract the third column of the velocity array.

Notes

This method is a functional alias for vz().

Returns

A 1D array of z-velocities in dimensionless internal units.

Return type

ndarray

cart2sph(pos=None)

Transform Cartesian coordinates (x, y, z) into spherical coordinates (r, phi, theta).

Convert 3D vectors into their spherical equivalents. The resulting coordinates follow the physics convention: - r: Radial distance. - p (phi): Azimuthal angle in the xy-plane [-pi, pi]. - t (theta): Polar angle (co-latitude) from the z-axis [0, pi].

Parameters

pos (ndarray, optional) – A (N, 3) array of Cartesian coordinates to transform. If None, use the instance’s self.pos attribute. Default is None.

Returns

A (N, 3) float32 array where columns represent [r, p, t] in dimensionless internal units.

Return type

ndarray

cm(ptype=None, r_max=None, center=array([0., 0., 0.]))

Calculate the center of mass of the model.

Parameters
  • ptype (int or list, optional) – Particle type(s) used to identify the center of mass minimum.

  • r_max (float, optional) – Maximum radius from the provided center to consider for the calculation.

  • center (array_like, optional) – The reference point (x, y, z) used if r_max is specified. Defaults to [0, 0, 0].

Returns

A 3D array [x, y, z] representing the center of mass coordinates.

Return type

ndarray

See also

get_potcenter

Get the potential minimum.

get_histocenter3D

Get the maximal density center using 3D histogram.

get_center_shrinking_sphere

Get the center using shrinking sphere.

cv

Get the velocity center of mass.

cv(ptype=None, r_max=None, center=array([0., 0., 0.]))

Calculate the velocity center of mass of the model.

Parameters
  • ptype (int or list, optional) – Particle type(s) used to identify the center of mass minimum.

  • r_max (float, optional) – Maximum radius, in positions space, from the provided center to consider for the calculation.

  • center (array_like, optional) – The reference point (x, y, z) used if r_max is specified. Defaults to [0, 0, 0].

Returns

A 3D array [x, y, z] representing the velocity center of mass coordinates.

Return type

ndarray

See also

cm

Get the center of mass.

get_potcenter

Get the potential minimum.

get_histocenter

Get the maximal density center using 1D histograms.

get_histocenter3D

Get the maximal density center using 3D histogram.

get_center_shrinking_sphere

Get the center using shrinking sphere.

dv_mean()

Estimate the mean relative velocity between particles.

Approximate the characteristic velocity difference between neighboring particles by scaling the total velocity dispersion by the particle count density.

Returns

Estimated mean relative speed between particles.

Return type

float

dx_mean()

Estimate the mean inter-particle separation.

Approximate the average distance between particles by treating the total number of bodies as being distributed within a characteristic volume defined by the system’s position dispersion.

Returns

Estimated mean physical distance between particles.

Return type

float

ekin()

Calculate the total specific kinetic energy of the system.

Computes the total kinetic energy per unit mass. This represents the energy-weighted velocity dispersion of the system.

Returns

Specific kinetic energy ($e_{kin}$) in internal units.

Return type

float

ellipticity_2D(axis='z', center=False)

Compute the 2D ellipticity of the particle distribution projection.

The ellipticity is calculated using the eigenvalues of the covariance matrix of the particle positions. This represents how “stretched” the distribution is in a 2D plane.

Parameters
  • axis ({'x', 'y', 'z'}, optional) – the axis to project along (i.e., the line of sight). ‘z’ projects onto the xy-plane, ‘y’ onto the xz-plane, and ‘x’ onto the yz-plane. Default is ‘z’.

  • center (bool, optional) – If True, subtracts the mean position (centroid) before computing the covariance. Default is False.

Returns

The ellipticity $epsilon = 1 - (b/a)$, where $a$ and $b$ are the major and minor semi-axes.

Return type

float

ellipticity_3D(center=False)

Compute the 3D shape parameters using the covariance matrix.

Quantifies the triaxiality of the particle distribution in 3D space based on the principal axes lengths ($a ge b ge c$).

Parameters

center (bool, optional) – If True, centers the distribution at the origin. Default is False.

Returns

  • ellipticity_1 (float) – Primary ellipticity defined as $1 - (b/a)$.

  • ellipticity_2 (float) – Secondary ellipticity defined as $1 - (c/a)$.

epot(eps)

Calculate the total specific potential energy.

Computes the total gravitational potential energy divided by the total mass of the system.

Parameters

eps (float) – The softening length (gravitational smoothing scale).

Returns

Specific potential energy ($w$ or $phi$) in internal units.

Return type

float

Warning

This method is currently serial and does not support MPI parallelization. It should only be called on a single process.

get_mass_tot()

Calculate the total mass of the system across all MPI processes.

Sum the individual masses of all particles in the current snapshot. This method uses a global MPI reduction to ensure the mass is aggregated across the entire distributed simulation.

Notes

This method is similar to totalmass().

Returns

  • float

  • | The total mass of the system in internal simulation units.

inertial_tensor()

Compute the rigid-body inertia tensor of the system.

The tensor is defined as: $I_{ij} = sum m (r^2 delta_{ij} - r_i r_j)$

This represents the distribution’s resistance to rotational acceleration about its principal axes.

Returns

A 3x3 symmetric matrix representing the inertia tensor.

Return type

ndarray

l()

Calculate the specific angular momentum vector for each particle.

Computes the angular momentum per unit mass, $mathbf{j} = mathbf{r} imes mathbf{v}$, for every particle.

Returns

A (N, 3) float array containing the specific angular momentum components $[j_x, j_y, j_z]$ for each particle.

Return type

ndarray

ltot()

Calculate the total specific angular momentum of the system.

Computes the mass-weighted average angular momentum of the system: $mathbf{j}_{tot} =

rac{sum m_i (mathbf{r}_i imes mathbf{v}_i)}{sum m_i}$.

ndarray

A (3,) float array representing the total specific angular momentum vector $[j_{x,tot}, j_{y,tot}, j_{z,tot}]$.

minert()

Calculate the root-mean-square (RMS) dispersion along Cartesian axes.

Computes the mass-weighted diagonal components of the second moment of the position distribution. This provides a measure of the “width” of the model along x, y, and z.

Returns

A 1D array [rms_x, rms_y, rms_z].

Return type

ndarray

phi_xy()

Calculate the cylindrical azimuthal angle.

Determine the polar angle in the xy-plane, ranging from the x-axis toward the y-axis.

Returns

A 1D array of azimuthal angles in radians in dimensionless units.

Return type

ndarray

phi_xyz()

Calculate the spherical azimuthal angle (longitude).

Determine the angle in the xy-plane relative to the x-axis, typically denoted as phi in spherical coordinates.

Returns

A 1D array of azimuthal angles in radians, ranging from -pi to pi.

Return type

ndarray

property r

Alias for rxyz().

rxy(center=None)

Calculate the 2D cylindrical radial distance.

Compute the projected distance from the origin or a specified center onto the xy-plane.

Parameters

center (array_like, optional) – A 3-element sequence (x, y, z) defining the origin. Only the x and y components are used. Default is None.

Returns

A 1D array of cylindrical radii in dimensionless internal units.

Return type

ndarray

rxyz(center=None)

Calculate the 3D spherical radial distance.

Compute the Euclidean distance from the origin or a specified center to each particle in 3D space.

Parameters

center (array_like, optional) – A 3-element sequence (x, y, z) defining the origin of the coordinate system. Default is None.

Returns

A 1D array of spherical radii in dimensionless internal units.

Return type

ndarray

size()

Estimate the characteristic size of the model using the moment of inertia.

Compute the spatial extent of the system based on the principal axes of the inertia tensor. Specifically, it returns the maximum value among the calculated moments of inertia, which serves as a robust proxy for the system’s overall radius or “spread.”

Returns

The maximum inertial scale of the particle distribution.

Return type

float

sph2cart(pos=None)

Transform spherical coordinates (r, phi, theta) into Cartesian coordinates (x, y, z).

Convert spherical components back into 3D Cartesian space. This method assumes the physics convention: - r: Radial distance. - p (phi): Azimuthal angle in the xy-plane. - t (theta): Polar angle (co-latitude) from the z-axis.

Parameters

pos (ndarray, optional) – A (N, 3) array of spherical coordinates [r, p, t]. If None, use the instance’s self.pos attribute (assuming it contains spherical data). Default is None.

Returns

A (N, 3) float32 array of Cartesian coordinates [x, y, z] in dimensionless internal units.

Return type

ndarray

theta_xyz()

Calculate the spherical elevation angle (latitude).

Compute the angle relative to the xy-plane (equator), typically denoted as theta in physics or latitude in geography.

Returns

A 1D array of elevation angles in radians, ranging from -pi/2 to pi/2.

Return type

ndarray

tork(acc)

Calculate the total torque acting on the system.

Computes the global torque vector $oldsymbol{ au} = sum m_i (mathbf{r}_i imes mathbf{a}_i)$ by aggregating the cross product of particle positions and their respective accelerations across all MPI processes.

Parameters

acc (ndarray) – A (N, 3) float array of accelerations $[a_x, a_y, a_z]$ for each particle.

Returns

A (3,) float array representing the total torque vector $[ au_x, au_y, au_z]$.

Return type

ndarray

Notes

This method uses a global MPI reduction to sum the torque contributions from all particles in the distributed simulation.

totalmass()

Calculate the total mass of the system across all MPI processes.

Sum the individual masses of all particles in the current snapshot. This method uses a global MPI reduction to ensure the mass is aggregated across the entire distributed simulation.

Notes

This method is similar to get_mass_tot().

Returns

  • float

  • | The total mass of the system in internal simulation units.

v_sigma()

Calculate the 3D root-mean-square (RMS) velocity dispersion.

Compute the mass-weighted standard deviation of the particle velocities relative to the center of velocity. This is a proxy for the kinetic temperature of the system.

Returns

The RMS velocity dispersion ($sigma_v$) in internal velocity units.

Return type

float

vn()

Calculate the magnitude (norm) of the velocity vectors.

Compute the Euclidean norm for each particle using the Cartesian velocity components.

Returns

A 1D array of velocity magnitudes in dimensionless internal units.

Return type

ndarray

vrxyz()

Calculate the radial velocity in the spherical coordinate system.

Compute the projection of the velocity vector onto the radial unit vector for each particle. This represents the rate of change of the spherical radius $r$.

Returns

A 1D array of radial velocities ($v_r$) in dimensionless internal units. Positive values indicate outward motion.

Return type

ndarray

vx()

Return the x-components of the velocity.

Extract the first column of the velocity array.

Returns

A 1D array of x-velocities in dimensionless internal units.

Return type

ndarray

vy()

Return the y-components of the velocity.

Extract the second column of the velocity array.

Returns

A 1D array of y-velocities in dimensionless internal units.

Return type

ndarray

vz()

Return the z-components of the velocity.

Extract the third column of the velocity array.

Returns

A 1D array of z-velocities in dimensionless internal units.

Return type

ndarray

x(center=None)

Return the x-coordinates of the particles.

Extract the first column of the position array. If a center is provided, offset the coordinates relative to that point.

Parameters

center (array_like, optional) – A 3-element sequence representing the (x, y, z) coordinates of the new origin. Default is None.

Returns

A 1D array of x-coordinates in dimensionless internal units.

Return type

ndarray

x_sigma()

Calculate the 3D root-mean-square (RMS) position dispersion.

Compute the mass-weighted standard deviation of the particle distances relative to the center of mass. This serves as a scalar measure of the spatial “spread” of the system.

Returns

The RMS position dispersion ($sigma_x$) in internal length units.

Return type

float

y(center=None)

Return the y-coordinates of the particles.

Extract the second column of the position array. If a center is provided, offset the coordinates relative to that point.

Parameters

center (array_like, optional) – A 3-element sequence representing the (x, y, z) coordinates of the new origin. Default is None.

Returns

A 1D array of y-coordinates in dimensionless internal units.

Return type

ndarray

z(center=None)

Return the z-coordinates of the particles.

Extract the third column of the position array. If a center is provided, offset the coordinates relative to that point.

Parameters

center (array_like, optional) – A 3-element sequence representing the (x, y, z) coordinates of the new origin. Default is None.

Returns

A 1D array of z-coordinates in dimensionless internal units.

Return type

ndarray