swiftsimio.snapshot_writer module

Provide utilities to create SWIFT files usable as initial conditions.

swiftsimio.snapshot_writer.get_dimensions(dimension: <module 'unyt.dimensions' from '/home/docs/checkouts/readthedocs.org/user_builds/swiftsimio/envs/stable/lib/python3.12/site-packages/unyt/dimensions.py'>) dict[source]

Return exponents corresponding to base dimensions for given unyt.dimensions.

Parameters:

dimension (unyt.dimensions) – Dimension for which we’re identifying the exponents.

Returns:

Dictionary of exponents corresponding to each base dimension.

Return type:

dict

Examples

>>> get_dimensions(unyt.dimensions.velocity)
{(mass): 0.0, (length): 1.0, (time): -1.0, (temperature): 0.0, (angle): 0.0,
 (current_mks): 0.0, 1: 0.0, (luminous_intensity): 0.0, (logarithmic): 0.0}
swiftsimio.snapshot_writer.generate_getter(name: str) Callable[[__SWIFTWriterParticleDataset], cosmo_array][source]

Generate a getter for the cosmo_array for name.

Parameters:

name (str) – Name of data field.

Returns:

Getter function that returns cosmo_array for name.

Return type:

Callable

swiftsimio.snapshot_writer.generate_setter(name: str, dimensions: <module 'unyt.dimensions' from '/home/docs/checkouts/readthedocs.org/user_builds/swiftsimio/envs/stable/lib/python3.12/site-packages/unyt/dimensions.py'>, unit_system: UnitSystem | str) Callable[[__SWIFTWriterParticleDataset, cosmo_array], None][source]

Generate a function that sets self._<name> to the value that is passed to it.

Parameters:
  • name (str) – String to set self._name to.

  • dimensions (unyt.dimensions) – Physical dimension of self._name (for consistency check).

  • unit_system (unyt.UnitSystem or str) – Unit system of self._name.

Returns:

Function to set self._name.

Return type:

Callable

swiftsimio.snapshot_writer.generate_deleter(name: str) Callable[[__SWIFTWriterParticleDataset], None][source]

Generate a function that destroys self._<name> (sets it back to None).

Parameters:

name (str) – Name of object to be destroyed.

Returns:

Function to delete self._<name>.

Return type:

Callable

swiftsimio.snapshot_writer.generate_dataset(writer: SWIFTSnapshotWriter, unit_system: UnitSystem | str, particle_type: int) __SWIFTWriterParticleDataset[source]

Generate a SWIFTWriterParticleDataset.

We _must_ do the following _outside_ of the class itself, as one can assign properties to a _class_ but not _within_ a class dynamically.

Here we loop through all of the possible properties in the metadata file. We then use the builtin property() function and some generators to create setters and getters for those properties. This will allow them to be accessed from outside by using SWIFTWriterParticleDataset(...).<name>, where the name is, for example, coordinates.

Parameters:
  • writer (SWIFTSnapshotWriter) – A reference to the writer that will contain the generated particle dataset writer.

  • unit_system (unyt.UnitSystem or str) – Unit system of the dataset.

  • particle_type (int) – The particle type of the dataset. Numbering convention is the same as SWIFT, with 0 corresponding to gas, etc. as usual.

Returns:

An empty dataset for the given particle type.

Return type:

SWIFTWriterParticleDataset

class swiftsimio.snapshot_writer.SWIFTSnapshotWriter(*, unit_system: UnitSystem | str = cosmological Unit System  Base Units:   length: Mpc   mass: 10000000000.0*Msun   time: 977.792221513146*Gyr   temperature: K   angle: rad   current_mks: A   luminous_intensity: cd   logarithmic: Np  Other Units:, boxsize: cosmo_array, dimension: int = 3, compress: bool = True, extra_header: dict | None = None, scale_factor: float32 = 1.0)[source]

Bases: object

The SWIFT dataset writer.

This is used to store all particle arrays and do some extra consistency checks and processing before writing a HDF5 file containing:

  • Fully consistent unit system

  • All required arrays for SWIFT to start

  • Required metadata (all automatic, apart from those required by __init__)

The written output can be read back in using swiftsimio, or used as initial conditions for a SWIFT simulation.

Any of the usual datasets (gas, dark_matter, stars, black_holes, sinks, neutrinos, boundary) can be populated with data. If a dataset is populated, it must have all required_fields filled in. The required fields can be viewed, for example, with swiftsimio.metadata.writer.required_fields.gas.keys(). The only exception is particle_ids, that will be generated automatically if missing. Smoothing length estimates for uniform-density gas in a periodic box can also be auto-generated, but this has to be done explicitly (see example below).

To fill in a required dataset, provide the data with a cosmo_array. For example:

w.gas.masses = cosmo_array(
    np.ones(n_p, dtype=float) * 1e6,
    u.solMass,
    comoving=True,
    scale_factor=w.scale_factor,
    scale_exponent=0,
)
Parameters:
  • unit_system (unyt.UnitSystem or str) – Unit system for dataset.

  • boxsize (cosmo_array) – Size of simulation box and associated units.

  • dimension (int, optional) – Dimensions (length of coordinate vector) of simulation.

  • compress (bool, optional) – Flag to turn on gzip compression.

  • extra_header (dict, optional) – Dictionary containing extra fields to write to the HDF5 Header group.

  • scale_factor (np.float32) – Scale factor associated with dataset.

Examples

An example showing populating all required fields for the gas dataset:

import numpy as np
import unyt as u
from swiftsimio import Writer, cosmo_array
from swiftsimio.metadata.writer.unit_systems import cosmo_unit_system

# number of gas particles
n_p = 1000
# scale factor of 1.0
a = 1.0
# Box is 100 Mpc
lbox = 100
boxsize = cosmo_array(
     [lbox, lbox, lbox],
     u.Mpc,
     comoving=True,
     scale_factor=a,
     scale_exponent=1,
)

# Create the Writer object. cosmo_unit_system corresponds to default Gadget-like
# units of 10^10 Msun, Mpc, and km/s
w = Writer(unit_system=cosmo_unit_system, boxsize=boxsize, scale_factor=a)

# Randomly spaced coordinates from 0 to lbox Mpc in each direction
w.gas.coordinates = cosmo_array(
    np.random.rand(n_p, 3) * lbox,
    u.Mpc,
    comoving=True,
    scale_factor=w.scale_factor,
    scale_exponent=1,
)

# Random velocities from 0 to 1 km/s
w.gas.velocities = cosmo_array(
    np.random.rand(n_p, 3),
    u.km / u.s,
    comoving=True,
    scale_factor=w.scale_factor,
    scale_exponent=1,
)

# Generate uniform masses as 10^6 solar masses for each particle
w.gas.masses = cosmo_array(
    np.ones(n_p, dtype=float) * 1e6,
    u.solMass,
    comoving=True,
    scale_factor=w.scale_factor,
    scale_exponent=0,
)

# Generate internal energy corresponding to 10^4 K
w.gas.internal_energy = cosmo_array(
    np.ones(n_p, dtype=float) * 1e4 / 1e6,
    u.kb * u.K / u.solMass,
    comoving=True,
    scale_factor=w.scale_factor,
    scale_exponent=-2,
)

# Generate initial guess for smoothing lengths based on mean inter-particle spacing
w.gas.generate_smoothing_lengths()

# w.gas.particle_ids can optionally be set, otherwise they are auto-generated

# write the initial conditions out to a file
w.write("ics.hdf5")
create_particle_datasets() None[source]

Create particle dataset for each particle type in the metadata.

write(filename: str) None[source]

Write the information in the dataset to file.

Parameters:

filename (str) – File to write to.