Methods description
- H5CosmoKit.Pk_suffix(ptype)
Maps a particle type to its corresponding label.
- Parameters:
ptype (list) – Particle type identifier(s), expected values are [0], [1], [4], [5], or [0,1,4,5].
- Returns:
Label of the particle type (‘g’ for gas, ‘c’ for cold dark matter, ‘s’ for stars, ‘bh’ for black holes, ‘m’ for all combined).
- Return type:
str
- Raises:
Exception – If no label is found for the provided ptype.
- H5CosmoKit.calc_entropy(U, rho_g)
Computes the entropy of gas particles.
- Parameters:
U (numpy.ndarray) – Specific internal energy array from snapshot (in km/s)^2.
rho_g (numpy.ndarray) – Density array from snapshot.
- Returns:
Computed entropy.
- Return type:
numpy.ndarray
- H5CosmoKit.calc_pressure(U, rho_g)
Computes the pressure of gas particles.
- Parameters:
U (numpy.ndarray) – Specific internal energy array from snapshot (in km/s)^2.
rho_g (numpy.ndarray) – Density array from snapshot.
- Returns:
Computed pressure.
- Return type:
numpy.ndarray
- H5CosmoKit.calc_soundSpeed(U)
Calculates the sound speed for an ideal gas given the specific internal energy.
The formula used assumes an adiabatic index (gamma) of 5/3. The sound speed is calculated using the simplified formula derived under the assumption of a constant volume (isochoric) process, where the density dependency can be ignored.
- Parameters:
U (numpy.ndarray) – An array of specific internal energies for which the sound speeds are calculated.
- Returns:
An array of sound speeds corresponding to the input specific internal energies.
- Return type:
numpy.ndarray
Note
This function is derived from the more general form: Cs = sqrt(gamma * P / rho) Where P = (gamma - 1) * rho * U is the pressure of the gas. In contexts where density (rho) is constant or its effect is normalized out, the formula simplifies to only depend on the specific internal energy (U): Cs = sqrt(gamma * (gamma - 1) * U)
- H5CosmoKit.calc_temperature(U, ne)
Computes the temperature of particles in a snapshot.
This function calculates the temperature of particles in a snapshot based on their internal energy, electron abundance, and the helium mass fraction. It uses constants for the Boltzmann constant and proton mass.
- Parameters:
U (numpy.ndarray) – Particle internal energy in (km/s)^2.
ne (numpy.ndarray) – Electron abundance.
- Returns:
An array of temperatures for each particle in the snapshot.
- Return type:
numpy.ndarray
- H5CosmoKit.compute_power_spectrum(pos, mass, BoxSize, grid=512, MAS='CIC')
Computes the power spectrum for given positions and masses of particles.
- H5CosmoKit.download_file(url, local_filename)
Downloads a file from a specified URL and saves it to a local path.
This function retrieves a file from the given URL and writes it to a local file, handling it in chunks to manage memory efficiently.
- Parameters:
url (str) – The URL of the file to download.
local_filename (str) – The local path (including filename) where the file will be saved.
- Returns:
The path to the downloaded file.
- Return type:
str
- Raises:
HTTPError – An error occurs from the HTTP request like 404, 500, etc.
- H5CosmoKit.interpolate_quantity(pos_g, quantity_g, boxSize)
Interpolates snapshot quantity onto a grid.
- Parameters:
pos_g (numpy.ndarray) – Particle coordinates in Mpch^{-1}.
quantity_g (numpy.ndarray) – Particle quantity to be interpolated.
boxSize (float) – Box size in Mpch^{-1}.
- Returns:
Interpolation function.
- Return type:
scipy.interpolate.NearestNDInterpolator
- H5CosmoKit.plot_internalenergy_distribution(path, snapshot_numbers, bw=7, x_limits=(0, 200000), sample_size=None, snapshot_base_name='snapshot')
Plots the distribution of internal energy from multiple snapshot files as a raincloud plot.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_base_name (str) – Base name of the snapshot files. Default is ‘snapshot’.
snapshot_numbers (list of int) – List of snapshot numbers to plot. Default is None.
bw (float) – Bandwidth for the density estimation in the raincloud plot. Default is 7.
x_limits (tuple) – X-axis limits for the plot. Default is (0, 100000).
sample_size (int, optional) – Number of random samples to select from each snapshot. If None, all data is used. Default is None.
Example
- plot_internalenergy_distribution(
path=’/gpfs/data/fs72085/mfo/CAMELS/CV0’, snapshot_numbers=[32, 44, 60, 90], sample_size=50000
)
This will display a raincloud plot of internal energies across various redshifts, with mean and median trends.
- H5CosmoKit.plot_median_internalenergy_with_polynomial_fit(path, snapshot_numbers, max_degree=5, snapshot_base_name='snapshot')
Plots the median internal energy against the scale factor from snapshot files and fits a piecewise polynomial.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_numbers (list of int) – List of snapshot numbers to be processed.
max_degree (int) – Maximum degree of the polynomial fit. Default is 5.
snapshot_base_name (str) – Base name of the snapshot files. Default is ‘snapshot’.
Example: plot_median_internalenergy_with_polynomial_fit(
path=’/gpfs/data/fs72085/mfo/CAMELS/CV0’, snapshot_numbers=[14, 18, 24, 28, 32, 34, 36, 38, 40, 42, 44, 46, 48, 50, 52, 54, 56, 58, 60, 62, 64, 66, 68, 70, 72, 74, 76, 78, 80, 82, 84, 86, 88, 90], max_degree=5
)
- H5CosmoKit.plot_median_soundspeed_with_polynomial_fit(path, snapshot_numbers, max_degree=5, a_0=0.1, snapshot_base_name='snapshot')
Plots the median sound speed against the scale factor from snapshot files and fits a polynomial.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_numbers (list of int) – List of snapshot numbers to be processed.
max_degree (int) – Maximum degree of the polynomial fit. Default is 5.
a_0 (float) – Scale factor threshold for the piecewise function.
snapshot_base_name (str) – Base name of the snapshot files. Default is ‘snapshot’.
Example: plot_median_soundspeed_with_polynomial_fit(
path=’/gpfs/data/fs72085/mfo/CAMELS/CV0’, snapshot_numbers=[14, 18, 24, 28, 32, 34, 36, 38, 40, 42, 44, 46, 48, 50, 52, 54, 56, 58, 60, 62, 64, 66, 68, 70, 72, 74, 76, 78, 80, 82, 84, 86, 88, 90], max_degree=5, a_0=0.1
)
- H5CosmoKit.plot_phase_diagram(ax, rho_g, quantity, title, xlabel, ylabel, bins=100, cmap='Reds')
Plots a 2D histogram phase diagram for given density and quantity (e.g., temperature, pressure, entropy).
- Parameters:
ax (matplotlib.axes.Axes) – The axis to plot on.
rho_g (numpy.ndarray) – Density array from snapshot.
quantity (numpy.ndarray) – Computed quantity (e.g., temperature, pressure, entropy).
title (str) – Plot title.
xlabel (str) – Label for the x-axis.
ylabel (str) – Label for the y-axis.
bins (int) – Number of bins for the 2D histogram.
cmap (str) – Colormap to use.
- H5CosmoKit.plot_power_spectrum_ratio(k_baryon, Pk_baryon, k_dm, Pk_dm)
Plots the ratio of power spectra.
- H5CosmoKit.plot_snapshot(ax, grid_quantity, boxSize, title, quantity)
Plots a snapshot quantity grid.
- Parameters:
ax (matplotlib.axes.Axes) – The matplotlib axis to plot on.
grid_quantity (numpy.ndarray) – Grid of snapshot quantity.
boxSize (float) – Box size in Mpch^{-1}.
title (str) – Title for the plot.
quantity (str) – Quantity name (e.g., ‘gas_density’ or ‘gas_temperature’).
- H5CosmoKit.plot_soundspeed_distribution(path, snapshot_numbers, bw=2, x_limits=(0, 300), sample_size=None, snapshot_base_name='snapshot')
Plots the distribution of sound speeds from multiple snapshot files as a raincloud plot.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_base_name (str) – Base name of the snapshot files. Default is ‘snapshot’.
snapshot_numbers (list of int) – List of snapshot numbers to plot. Default is None.
bw (float) – Bandwidth for the density estimation in the raincloud plot. Default is 1.
x_limits (tuple) – X-axis limits for the plot. Default is (0, 300).
sample_size (int, optional) – Number of random samples to select from each snapshot. If None, all data is used. Default is None.
Example
- plot_soundspeed_distribution(
path=’/gpfs/data/fs72085/mfo/CAMELS/CV0’, snapshot_numbers=[32, 44, 60, 90], sample_size=50000
)
This will display a raincloud plot of sound speeds across various redshifts, with mean and median trends.
- H5CosmoKit.power_ratio(f_snap)
Processes a snapshot file to compute and plot power spectrum ratios.
- H5CosmoKit.preview(path, snapshot_numbers, quantity, snapshot_base_name='snapshot', unit_scale='kpc')
Creates a plot of snapshot quantity from multiple snapshot files.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_numbers (list of int) – List of snapshot numbers to plot.
quantity (str) – Quantity name (e.g., ‘gas_density’ or ‘gas_temperature’).
snapshot_base_name (str) – Base name of the snapshot files (default: ‘snapshot’). Underscore is accounted for.
unit_scale (str) – Unit scale (‘kpc’ or ‘mpc’) for interpreting the data. Default is ‘kpc’.
Example
from H5CosmoKit import preview
# Example usage to visualize gas density from a snapshot in a cosmological simulation path = ‘/your/data/directory’ snapshot_numbers = [10, 20, 30] # Example with multiple snapshots quantity = ‘gas_density’
preview(path, snapshot_numbers, quantity, unit_scale=’Mpc’)
- H5CosmoKit.preview_3d(path, snapshot_numbers, quantity, subset_size, snapshot_base_name='snapshot', unit_scale='kpc')
Creates an interactive 3D plot of snapshot quantity using Plotly and saves it as an HTML file.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_numbers (list of int) – List of snapshot numbers to plot.
quantity (str) – Quantity name (e.g., ‘gas_density’ or ‘gas_temperature’).
subset_size (int) – Size of randomly chosen points to be plotted to manage performance.
snapshot_base_name (str) – Base name of the snapshot files (default: ‘snapshot’). Underscore is accounted for.
Example
path = ‘/your/data/directory’ snapshot_numbers = [30] quantity = ‘gas_temperature’ subset_size = 10000
preview_3d(path, snapshot_numbers, quantity, subset_size, unit_scale=’Mpc’)
- H5CosmoKit.preview_phase_diagram(path, snapshot_numbers, quantity, snapshot_base_name='snapshot', unit_scale='kpc')
Generates phase diagrams for density vs temperature, pressure, or entropy from multiple snapshots.
- Parameters:
path (str) – The base path containing the snapshot files.
snapshot_numbers (list of int) – List of snapshot numbers to plot.
quantity (str) – Quantity to plot (‘temperature’, ‘pressure’, ‘entropy’).
snapshot_base_name (str) – Base name of the snapshot files.
unit_scale (str) – Unit scale (‘kpc’ or ‘mpc’) for interpreting the data.
- H5CosmoKit.read_particles(data, part_type, mass_table)
Reads positions and masses of particles of a given type from the HDF5 data.
- H5CosmoKit.read_snapshot(snapshot_path, unit_scale='kpc')
Reads snapshot data from a HDF5 file, converting units if required.
- Parameters:
snapshot_path (str) – Path to the HDF5 snapshot file.
unit_scale (str) – Unit scale for length (‘kpc’ or ‘mpc’) of your snapshot file, default is ‘mpc’.
- Returns:
- Contains the box size, redshift, scale_factor, positions, densities, internal energies,
and electron abundances (if available) adjusted to the desired unit scale.
- Return type:
tuple
System of units for default option (mpc): %—- System of units UnitLength_in_cm 3.085678e24 % 1.0 Mpc UnitMass_in_g 1.989e43 % 1.0e10 solar masses UnitVelocity_in_cm_per_s 1e5 % 1 km/sec
- H5CosmoKit.save_plot(fig, path, quantity)
Saves the figure to a file with an appropriate filename.
- Parameters:
fig (matplotlib.figure.Figure) – The figure to be saved.
path (str) – The base path where the plot will be saved.
quantity (str) – Quantity name (e.g., ‘gas_density’ or ‘gas_temperature’).
- H5CosmoKit.validate_quantity(quantity)
Validates the provided quantity against allowed values.
- Parameters:
quantity (str) – The quantity to validate.
- Returns:
The validated quantity.
- Return type:
str