swiftgalaxy package

Extend SWIFTSimIO to analyse individual galaxies.

SWIFTGalaxy is a module that extends SWIFTSimIO tailored to analyses of particles belonging to individual simulated galaxies. It inherits from and extends the functionality of the SWIFTDataset. It understands the content of halo catalogues (supported: SOAP, Velociraptor, Caesar) and therefore which particles belong to a galaxy or other group of particles, and its integrated properties. The particles occupy a coordinate frame that is enforced to be consistent, such that particles loaded on-the-fly will match e.g. rotations and translations of particles already in memory. Intuitive masking of particle datasets is also enabled. Finally, some utilities to make working in cylindrical and spherical coordinate systems more convenient are also provided.

class swiftgalaxy.Caesar(caesar_file: str | None = None, group_type: str | None = None, group_index: int | Sequence[int] | None = None, centre_type: str = 'minpot', extra_mask: str | MaskCollection = 'bound_only', custom_spatial_offsets: cosmo_array | None = None)[source]

Bases: _HaloCatalogue

Interface to Caesar halo catalogues for use with swiftgalaxy.

Takes a caesar output file and configuration options and provides an interface that swiftgalaxy understands. Also exposes the halo/galaxy properties computed by CAESAR for a single object of interest with the same interface provided by the Group class in the caesar python package. Reading of properties is done on-the-fly, and only rows corresponding to the object of interest are read from disk.

Parameters:
  • caesar_file (str) – The catalogue file (hdf5 format) output by caesar.

  • group_type (str) – The category of the object of interest, either "halo" or "galaxy".

  • group_index (int) – Position(s) of the object(s) of interest in the catalogue arrays. In the case of multiple objects, duplicate entries are not allowed.

  • centre_type (str (optional), default: "minpot") – Type of centre, chosen from those provided by caesar. Default is the position of the particle with the minimum potential, "minpot", alternatively "" can be used for the centre of mass.

  • extra_mask (str or MaskCollection (optional), default: "bound_only") – Mask to apply to particles after spatial masking. If "bound_only", then the galaxy is masked to include only the gravitationally bound particles as provided by caesar. A user-defined mask can also be as an an object (such as a swiftgalaxy.masks.MaskCollection) that has attributes with names corresponding to present particle names (e.g. gas, dark_matter, etc.), each containing a mask.

  • custom_spatial_offsets (cosmo_array (optional), default: None) – A region to override the automatically-determined region enclosing group member particles. May be used in conjunction with extra_mask, for example to select all simulation particles in an aperture around the object of interest (see ‘Masking’ section of documentation for a cookbook example). Provide a pair of offsets from the object’s centre along each axis to define the region, for example for a cube extending +/- 1 Mpc from the centre: cosmo_array([[-1.0, 1.0], [-1.0, 1.0], [-1.0, 1.0]], u.Mpc).

Notes

Note

CAESAR only supports index access to catalogue lists, not identifier access. This means that the group_index is simply the position of the object of interest in the catalogue list.

Examples

Given a file s12.5n128_0012.hdf5 at /output/path/, the following creates a Caesar object for the entry at index 3 in the catalogue (i.e. the 4th row, indexed from 0) and demonstrates retrieving its virial mass.

>>> cat = Caesar(
>>>     caesar_file="/output/path/s12.5n128_0012.hdf5",
>>>     group_type="halo",
>>>     group_index=3,  # 4th entry in catalogue (indexed from 0)
>>> )
>>> cat.info()
{'GroupID': 3,
'ages': {'mass_weighted': unyt_quantity(2.26558173, 'Gyr'),
         'metal_weighted': unyt_quantity(2.21677032, 'Gyr')},
'bh_fedd': unyt_quantity(3.97765937, 'dimensionless'),
'bhlist_end': 12,
'bhlist_start': 11,
...
'virial_quantities': {'circular_velocity': unyt_quantity(158.330253, 'km/s'),
                      'm200c': unyt_quantity(1.46414384e+12, 'Msun'),
                      'm2500c': unyt_quantity(8.72801239e+11, 'Msun'),
                      'm500c': unyt_quantity(1.23571772e+12, 'Msun'),
                      'r200': unyt_quantity(425.10320408, 'kpccm'),
                      'r200c': unyt_quantity(327.46600342, 'kpccm'),
                      'r2500c': unyt_quantity(118.77589417, 'kpccm'),
                      'r500c': unyt_quantity(228.03265381, 'kpccm'),
                      'spin_param': unyt_quantity(2.28429179,'s/(Msun*km*kpccm)'),
                      'temperature': unyt_quantity(902464.88453405, 'K')}}
>>> cat.virial_quantities
{'circular_velocity': unyt_quantity(158.330253, 'km/s'),
 'm200c': unyt_quantity(1.46414384e+12, 'Msun'),
 'm2500c': unyt_quantity(8.72801239e+11, 'Msun'),
 'm500c': unyt_quantity(1.23571772e+12, 'Msun'),
 'r200': unyt_quantity(425.10320408, 'kpccm'),
 'r200c': unyt_quantity(327.46600342, 'kpccm'),
 'r2500c': unyt_quantity(118.77589417, 'kpccm'),
 'r500c': unyt_quantity(228.03265381, 'kpccm'),
 'spin_param': unyt_quantity(2.28429179, 's/(Msun*km*kpccm)'),
 'temperature': unyt_quantity(902464.88453405, 'K')}
>>> cat.virial_quantities["m200c"]
unyt_quantity(1.46414384e+12, 'Msun')
caesar_file: str | None
property centre: cosmo_array

Obtain the centre specified by the centre_type from the halo catalogue.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The centre(s) of the object(s) of interest.

Return type:

cosmo_array

centre_type: str
property count: int

Number of galaxies in the catalogue.

When in multi-galaxy mode and not masked to a single row this can be >1.

Returns:

The number of galaxies in the catalogue.

Return type:

int

property group_index: List[int] | int

Get the index of the object of interest in the halo catalogue.

In multi-galaxy mode when no mask is active this gets the list of indices. This is just the position in the arrays stored in the Caesar halo or galaxy lists (according to the group_type given at initialization).

Returns:

The index or indices of the object(s) of interest in the halo catalogue.

Return type:

int or list

group_type: str
property velocity_centre: cosmo_array

Obtain the velocity centre from the halo catalogue.

It is specified by the velocity_centre_type.

Returns:

The velocity centre provided by the halo catalogue.

Return type:

cosmo_array

velocity_centre_type: str
class swiftgalaxy.MaskCollection(**kwargs: slice | EllipsisType | ndarray[tuple[Any, ...], dtype[_ScalarT]] | None | LazyMask)[source]

Bases: object

Barebones container for mask objects.

Takes a set of kwargs at initialisation and assigns their values to attributes of the object. Attempts to access a non-existent attribute returns None instead of raising an AttributeError.

This is intended to hold masks that can be applied to cosmo_array objects under the names of particle types (e.g. gas, dark_matter, etc.), but this is not checked or enforced.

Parameters:

**kwargs – Any items passed as kwargs will have their values passed to correspondingly named attributes of this object.

Notes

Note

The velociraptor.swift.swift module makes some use of a namedtuple called MaskCollection. These objects are not valid where swiftgalaxy functions expect a MaskCollection because namedtuple objects are immutable.

Examples

n_dm = 123  # suppose this is number of dark matter particles
# all these masks select all particles:
MaskCollection(
    gas=np.s_[...],
    dark_matter=np.ones(n_dm, dtype=bool),
    stars=None
)
combined_with(other_mask_collection: MaskCollection, *, sg: SWIFTGalaxy) MaskCollection[source]

Combine this MaskCollection with another.

data[this_mask.<type>.mask][other_mask.<type>.mask] and data[combined_mask.<type>.mask] are equivalent, where combined_mask = this_mask.combined_with(other_mask).

Parameters:
class swiftgalaxy.SOAP(soap_file: str | None = None, soap_index: int | Sequence[int] | None = None, extra_mask: str | MaskCollection = 'bound_only', centre_type: str = 'input_halos.halo_centre', velocity_centre_type: str = 'bound_subhalo.centre_of_mass_velocity', custom_spatial_offsets: cosmo_array | None = None)[source]

Bases: _HaloCatalogue

Interface to SOAP halo catalogues for use with swiftgalaxy.

Takes a set of SOAP output files and configuration options and provides an interface that swiftgalaxy understands. Also exposes the galaxy properties computed by SOAP for a single object of interest through the swiftsimio interface.

Parameters:
  • soap_file (str, default: None) – The filename of a SOAP catalogue file, possibly including the path.

  • soap_index (int or list, default: None) – The position (row) in the SOAP catalogue corresponding to the object of interest. Duplicate entries are not allowed.

  • extra_mask (str or MaskCollection (optional), default: "bound_only") – Mask to apply to particles after spatial masking. If "bound_only", then the galaxy is masked to include only the gravitationally bound particles as determined by SOAP. A user-defined mask can also be provided as an an object (such as a swiftgalaxy.masks.MaskCollection) that has attributes with names corresponding to present particle names (e.g. gas, dark_matter, etc.), each containing a mask.

  • centre_type (str (optional), default: "input_halos.halo_centre") – Type of centre, chosen from those provided by SOAP. This should be expressed as a string analogous to what would be written in swiftsimio code (or swiftgalaxy) to access that property in the SOAP catalogue. The default takes the "input_halos.halo_centre" (usually the centre of potential, e.g. the HBT+ halo finder defines it in this way; another option amongst many more is "bound_subhalo.centre_of_mass".

  • velocity_centre_type (str (optional), default: "bound_subhalo.centre_of_mass_velocity") – Type of velocity centre, chosen from those provided by SOAP. This should be expressed as a string analogous to what would be written in swiftsimio code (or swiftgalaxy) to access that property in the SOAP catalogue. The default takes the "bound_subhalo.centre_of_mass_velocity"; note that there is no velocity corresponding to the centre of potential ("input_halos.halo_centre_velocity" is not defined). Another useful option could be "exclusive_sphere_1kpc.centre_of_mass_velocity" to choose the velocity of bound particles in the central 1 kpc.

  • custom_spatial_offsets – A region to override the automatically-determined region enclosing group member particles. May be used in conjunction with extra_mask, for example to select all simulation particles in an aperture around the object of interest (see ‘Masking’ section of documentation for a cookbook example). Provide a pair of offsets from the object’s centre along each axis to define the region, for example for a cube extending +/- 1 Mpc from the centre: cosmo_array([[-1.0, 1.0], [-1.0, 1.0], [-1.0, 1.0]], u.Mpc).

Notes

Note

SOAP only supports index access to catalogue arrays, not identifier access. This means that the soap_index is simply the position of the object of interest in the SOAP catalogue arrays.

Examples

Given a file halo_properties_0050.hdf5, the following creates a SOAP object for the entry at index 0 in the catalogue (i.e. the 1st row, indexed from 0) and demonstrates retrieving its virial mass (M200crit).

>>> cat = SOAP(
>>>     soap_file="/output/path/halo_properties_0050.hdf5"
>>>     soap_index=0,  # 1st entry in catalogue (indexed from 0)
>>> )
>>> cat.spherical_overdensity_200_crit.total_mass.to(u.Msun)
cosmo_array([6.72e+12], dtype=float32, units='1.98841586e+30*kg', comoving=False)
bound_subhalo: __SWIFTGroupDataset
property centre: cosmo_array

Obtain the centre from the halo catalogue.

It is specified by the centre_type.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The centre(s) of the object(s) of interest.

Return type:

cosmo_array

centre_type: str
property count: int

Number of galaxies in the catalogue.

When in multi-galaxy mode and not masked to a single row this can be >1.

Returns:

The number of galaxies in the catalogue.

Return type:

int

input_halos: __SWIFTGroupDataset
soap_file: str
property soap_index: int | List[int]

Get the index of the objects of interest in the halo catalogue.

In multi-galaxy mode when no mask is active this gets the list of indices.

This is just the position in the arrays stored in the SOAP SWIFTDataset (that could have a SWIFTMask).

Returns:

The index or indices of the object(s) of interest in the halo catalogue.

Return type:

int or list

property velocity_centre: cosmo_array

Obtain the velocity centre from the halo catalogue.

It is specified by the velocity_centre_type.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The centre(s) of the object(s) of interest.

Return type:

cosmo_array

velocity_centre_type: str
class swiftgalaxy.SWIFTGalaxies(snapshot_filename: str, halo_catalogue: _HaloCatalogue, auto_recentre: bool = True, preload: Set[str] = {}, transforms_like_coordinates: Set[str] = {}, transforms_like_velocities: Set[str] = {}, id_particle_dataset_name: str = 'particle_ids', coordinates_dataset_name: str = 'coordinates', velocities_dataset_name: str = 'velocities', coordinate_frame_from: SWIFTGalaxy | None = None, optimize_iteration: str = 'auto')[source]

Bases: object

Facilitates efficiently iterating over many objects of interest from a simulation.

SWIFT simulation snapshots contain particles grouped by “top-level cells” that cover the simulation volume. The minimum number of particles that it makes sense to read is therefore those contained in one such top-level cell. If one wants to create many SWIFTGalaxy objects from one simulation snapshot, there is a risk that the same data are read many times, such as when multiple target objects lie within the same top-level cell. This class provides a convenient way to iterate over multiple target objects while minimizing the I/O overhead by managing the order of iteration to group together target objects that occupy common top-level cells and only reading the data once.

An important consequence to be aware of is that the iteration order is not controlled by the user because it must be chosen to group objects in the same top-level cell(s) together. The iteration order is available as the iteration_order attribute of a SWIFTGalaxies object. Alternatively, output of a function applied to a list of target objects in the same order as the input list can be obtained using the map() method.

There is an obvious opportunity to parallelize the iteration process by passing each region (potentially each containing multiple target objects) to worker processes as they become available, for example. This current initial version of the SWIFTGalaxies class does not yet support parallel iteration, instead prioritizing the release of a working serial implementation. Support for parallelization will be added later as a high priority.

Parameters:
  • snapshot_filename (str) – Name of file containing snapshot.

  • halo_catalogue (_HaloCatalogue) – A halo catalogue instance from swiftgalaxy.halo_catalogues, e.g. a swiftgalaxy.halo_catalogues.SOAP instance. It should specify more than one target object, e.g. by setting its soap_index=[0, 123, 456, ...].

  • auto_recentre (bool (optional), default: True) – If True, the coordinate system will be automatically recentred on the position and velocity centres defined by the halo_catalogue.

  • preload (set (optional), default: set()) – Deprecated and ignored.

  • transforms_like_coordinates (set (optional), default: set()) – Names of fields that behave as velocities. It is assumed that these exist for all present particle types. When the coordinate system is rotated or boosted, the associated arrays will be transformed accordingly. The velocities dataset (or its alternative name given in the velocities_dataset_name parameter) is implicitly assumed to behave as velocities.

  • transforms_like_velocities (set (optional), default: set()) – Names of fields that behave as velocities. It is assumed that these exist for all present particle types. When the coordinate system is rotated or boosted, the associated arrays will be transformed accordingly. The velocities dataset (or its alternative name given in the velocities_dataset_name parameter) is implicitly assumed to behave as velocities.

  • id_particle_dataset_name (str (optional), default: "particle_ids") – Name of the dataset containing the particle IDs, assumed to be the same for all present particle types.

  • coordinates_dataset_name (str (optional), default: "coordinates") – Name of the dataset containing the particle spatial coordinates, assumed to be the same for all present particle types.

  • velocities_dataset_name (str (optional), default: "velocities") – Name of the dataset containing the particle velocities, assumed to be the same for all present particle types.

  • coordinate_frame_from (SWIFTGalaxy (optional), default: None) – Another SWIFTGalaxy to copy the coordinate frame (centre and rotation) and velocity coordinate frame (boost and rotation) from.

  • optimize_iteration (str (optional), default: "auto") – Can be "auto", "dense" or "sparse". See docstrings of methods _eval_sparse_optimized_solution() and _eval_dense_optimized_solution() for explanations of optimization schemes. In most cases leave set to default "auto" to automatically determine optimal solution.

Examples

Using SWIFTGalaxies is almost the same as using the main SWIFTGalaxy class, except that (i) the halo catalogue is initialized with multiple target objects and (ii) the SWIFTGalaxies class provides an iteration method (__iter__), and determines its own iteration order. For example:

from swiftgalaxy import SWIFTGalaxies, SOAP
sgs = SWIFTGalaxies(
    "snapshot.hdf5",
    SOAP(
        "soap.hdf5",
        soap_index=[0, 123, 456],  # multiple target indices
    ),
)
iteration_order = sgs.iteration_order  # be aware of the order of iteration
for sg in sgs:
    # some analysis involving the pre-loaded data fields goes here:
    sg.element_abundances.carbon
    sg.dark_matter.coordinates
    sg.stars.velocities

Alternatively the map() method can be used to apply a function to all of the SWIFTGalaxy’s created by this class. For example:

from swiftgalaxy import SWIFTGalaxies, SOAP
sgs = SWIFTGalaxies(
    "snapshot.hdf5",
    SOAP(
        "soap.hdf5",
        soap_index=[0, 123, 456],  # multiple target indices
    ),
)

def analysis(sg):
    # this function can also have additional args & kwargs, if needed
    # it should only access the pre-loaded data fields
    sg.element_abundances.carbon
    sg.dark_matter.coordinates
    sg.stars.velocities
    return sg.element_abundances.carbon.mean()

# map accepts arguments `args` and `kwargs`, passed through to function, if needed
result = sgs.map(analysis)
count: int
halo_catalogue: _HaloCatalogue
property iteration_order: ndarray

Property holding the order that the target objects will be iterated in.

The iteration order is likely not the same as the order that the targets are provided in because this is probably not an optimal iteration order. This property attribute provides the optimized iteration order evaluated by SWIFTGalaxies.

Returns:

Array of indices specifying the iteration order.

Return type:

numpy.ndarray

map(func: Callable, args: List[Tuple] | None = None, kwargs: List[Dict] | None = None) List[Any][source]

Apply a function to each object of interest and return a list of results.

The iteration order of SWIFTGalaxies is not necessarily the order that the objects of interest are provided by the user because the class determined an efficient iteration order to minimize I/O operations. This method applies a provided function to each object of interest in an efficient order then returns the results in a list ordered in the same order that the objects of interest were input.

The function to be evaluated should expect a SWIFTGalaxy (from those to be iterated over) as its first argument. It may accept lists of additional arguments and/or keyword arguments (with each element corresponding to one entry in the list of target objects) that can be passed to map as a tuple of arguments and a dict of keyword arguments.

Currently this function only executes serially but adding a parallel execution option, and further support for parallelization in analysis, is a high priority.

Parameters:
  • func (callable) – The function to be evaluated.

  • args (list (optional), default: None) – List of additional arguments to the function to be evaluated (the first argument is always the current SWIFTGalaxy in the iteration). Each item in the list should be a tuple of arguments, with one tuple for each galaxy being iterated over. See examples section for further details.

  • kwargs (list (optional), default: None) – List of additional keyword arguments to pass to the function to be evaluated. Each item in the list should be a dict of keyword arguments, with one dict for each galaxy being iterated over. Dictionary keys are the names of the keyword arguments and the corresponding dictionary values are the values of the keyword arguments. See examples section for further details.

Returns:

A list containing the return value(s) of the function applied to each object of interest, in the same order as the objects of interest were passed to the halo finder interface.

Return type:

list

Examples

A simple example that applies a function dm_median_position to each galaxy in a list of targets [11, 22, 33]:

from swiftgalaxy import SWIFTGalaxies, SOAP

# define the function that we will apply to each SWIFTGalaxy object:
def dm_median_position(sg):
    return np.median(sg.dark_matter.coordinates, axis=0)

sgs = SWIFTGalaxies(
    "my_snapshot.hdf5",
    SOAP(
        "my_soap.hdf5",
        soap_index=[11, 22, 33],
    ),
)
my_result = sgs.map(dm_median_position)

The result stored in my_result contains the result of the function for the galaxies at index 11, 22 and 33, in the same order as they are given in the soap_index list.

This second example shows how to pass extra arguments and/or keyword arguments to the function given to map:

from swiftgalaxy import SWIFTGalaxies, SOAP

# define the function that we will apply to each SWIFTGalaxy object:
def dm_median_position(
    sg,  # the first argument is always a SWIFTGalaxy from the iteration
    extra_argument_1,
    extra_argument_2,
    extra_kwarg_1=None,
    extra_kwarg_2=None,
):
    # presumably make use of the extra arguments and/or kwargs here...
    return np.median(sg.dark_matter.coordinates, axis=0)

sgs = SWIFTGalaxies(
    "my_snapshot.hdf5",
    SOAP("my_soap.hdf5",
    soap_index=[11, 22, 33]),
)
my_result = sg.map(
    dm_median_position,
    args=[
        (my_extra_arg_1_for_galaxy_11, my_extra_arg_2_for_galaxy_11),
        (my_extra_arg_1_for_galaxy_22, my_extra_arg_2_for_galaxy_22),
        (my_extra_arg_1_for_galaxy_33, my_extra_arg_2_for_galaxy_33),
    ],
    kwargs=[
        dict(
            extra_kwarg_1=my_extra_kwarg_1_for_galaxy_11,
            extra_kwarg_2=my_extra_kwarg_2_for_galaxy_11,
        ),
        dict(
            extra_kwarg_1=my_extra_kwarg_1_for_galaxy_22,
            extra_kwarg_2=my_extra_kwarg_2_for_galaxy_22,
        ),
        dict(
            extra_kwarg_1=my_extra_kwarg_1_for_galaxy_33,
            extra_kwarg_2=my_extra_kwarg_2_for_galaxy_33,
        ),
    ]
)

Note that if you have only a single extra argument it must still be packaged as a tuple, for instance:

args=[
    (my_extra_arg_for_galaxy_11, ),
    (my_extra_arg_for_galaxy_22, ),
    (my_extra_arg_for_galaxy_33, ),
]

The commas inside the parentheses are not optional!

class swiftgalaxy.SWIFTGalaxy(snapshot_filename: str, halo_catalogue: _HaloCatalogue | None, auto_recentre: bool = True, transforms_like_coordinates: Set[str] = {}, transforms_like_velocities: Set[str] = {}, id_particle_dataset_name: str = 'particle_ids', coordinates_dataset_name: str = 'coordinates', velocities_dataset_name: str = 'velocities', coordinate_frame_from: SWIFTGalaxy | None = None, _data_server: SWIFTGalaxy | None = None)[source]

Bases: SWIFTDataset

A representation of a simulated galaxy.

A SWIFTGalaxy represents a galaxy from a simulation, including both its particles and integrated properties. A halo finder catalogue is required to define which particles belong to the galaxy and to provide integrated properties. The implementation is an extension of the SWIFTDataset class, so all the functionality of such a dataset is also available for a SWIFTGalaxy. The swiftsimio.reader.__SWIFTGroupDataset objects familiar to swiftsimio users (e.g. a GasDataset) have an analogous _SWIFTGroupDatasetHelper class (e.g. GasDatasetHelper) that maintains their usual functionality and extends it with new features. swiftsimio.reader.__SWIFTNamedColumnDataset instances are have analogues as _SWIFTNamedColumnDatasetHelper objects.

For an overview of available features see the examples below, and the narrative documentation pages.

Parameters:
  • snapshot_filename (str) – Name of file containing snapshot.

  • halo_catalogue (_HaloCatalogue (optional), default: None) – A halo catalogue instance from swiftgalaxy.halo_catalogues, e.g. a swiftgalaxy.halo_catalogues.SOAP instance.

  • auto_recentre (bool (optional), default: True) – If True, the coordinate system will be automatically recentred on the position and velocity centres defined by the halo_catalogue.

  • transforms_like_coordinates (set (optional), default: set()) – Names of fields that behave as spatial coordinates. It is assumed that these exist for all present particle types. When the coordinate system is rotated or translated, the associated arrays will be transformed accordingly. The coordinates dataset (or its alternative name given in the coordinates_dataset_name parameter) is implicitly assumed to behave as spatial coordinates.

  • transforms_like_velocities (set (optional), default: set()) – Names of fields that behave as velocities. It is assumed that these exist for all present particle types. When the coordinate system is rotated or boosted, the associated arrays will be transformed accordingly. The velocities dataset (or its alternative name given in the velocities_dataset_name parameter) is implicitly assumed to behave as velocities.

  • id_particle_dataset_name (str (optional), default: "particle_ids") – Name of the dataset containing the particle IDs, assumed to be the same for all present particle types.

  • coordinates_dataset_name (str (optional), default: "velocities") – Name of the dataset containing the particle spatial coordinates, assumed to be the same for all present particle types.

  • velocities_dataset_name (str (optional), default: "velocities") – Name of the dataset containing the particle velocities, assumed to be the same for all present particle types.

  • coordinate_frame_from (SWIFTGalaxy (optional), default: None) – Another SWIFTGalaxy to copy the coordinate frame (centre and rotation) and velocity coordinate frame (boost and rotation) from.

  • _data_server (SWIFTGalaxy) – SWIFTGalaxy to use for data read operations. For example SWIFTGalaxies uses this to share read data between grouped objects.

See also

_SWIFTGroupDatasetHelper, _SWIFTNamedColumnDatasetHelper, swiftgalaxy.halo_catalogues

Examples

Assuming we have a snapshot file snap.hdf5, and velociraptor outputs halos.properties, halos.catalog_groups, etc., with the default names for coordinates, velocities and particle_ids, we can initialise a SWIFTGalaxy for the first row (indexed from 0) in the halo catalogue very easily:

from swiftgalaxy import SWIFTGalaxy, Velociraptor
mygalaxy = SWIFTGalaxy(
    'snap.hdf5',
    Velociraptor(
        'halos',
        halo_index=0
    )
)

Like a SWIFTDataset, the particle datasets are accessed as below, and all data are loaded ‘lazily’, on demand.

mygalaxy.gas.particle_ids
mygalaxy.dark_matter.coordinates

However, information from the halo catalogue is used to select only the particles identified as bound to this galaxy. The coordinate system is centred in both position and velocity on the centre and peculiar velocity of the galaxy, as determined by the halo catalogue. The coordinate system can be further manipulated, and all particle arrays will stay in a consistent reference frame at all times.

Again like for a SWIFTDataset, the units and metadata are available:

mygalaxy.units
mygalaxy.metadata

The halo catalogue interface is accessible as shown below. What this interface looks like depends on the halo catalogue being used, but will provide values for the individual galaxy of interest.

mygalaxy.halo_catalogue

In this case with Velociraptor, we can get the virial mass like this:

mygalaxy.halo_catalogue.masses.mvir

For a complete description of available features see the narrative documentation pages.

boost(boost: cosmo_array) None[source]

Apply a ‘boost’ to the velocity coordinates.

The provided velocity vector is added to all particle velocities. If the ‘reference velocity’ that should be set to zero is known, use recentre_velocity() instead. All datasets specified in the transforms_like_velocities argument to SWIFTGalaxy are transformed (by default velocities for all present particle types).

Parameters:

boost (cosmo_array) – The velocity to boost by.

property centre: cosmo_array

Get the current origin of the coordinate reference frame.

It is given with respect to the native simulation coordinate reference frame.

Returns:

The origin of the coordinate reference frame.

Return type:

cosmo_array

coordinates_dataset_name: str
create_datasets() None

Create datasets for present groups.

Present groups are specified in metadata.present_group_names.

These can then be accessed using their underscore names, e.g. gas.

filename: Path
get_metadata() None

Load the metadata from the SWIFT snapshot.

Ordinarily this happens automatically, but you can call this function again if you mess things up.

get_units() None

Load the units from the SWIFT snapshot.

Ordinarily this happens automatically, but you can call this function again if you mess things up.

halo_catalogue: _HaloCatalogue | None
property handle: File

Provide a property that returns the handle, opening it if necessary.

Returns:

The file handle.

Return type:

h5py.File

Raises:

RuntimeError – If the handle is not managed by this object but is falsy (e.g. it is a closed handle).

id_particle_dataset_name: str
mask_particles(mask_collection: MaskCollection) None[source]

Select a subset of the currently selected particles.

The masks to be applied can by in any format accepted by a cosmo_array via __getitem__() and should be collected in a swiftgalaxy.masks.MaskCollection. The selection is applied permanently to all particle datasets for this galaxy. Temporary masks (e.g. for interactive use) can be created by using the __getitem__() (square brackets) method of the SWIFTGalaxy, any of its associated _SWIFTGroupDatasetHelper or _SWIFTNamedColumnDatasetHelper objects, but note that to ensure internal consistency, these return a masked copy of the entire SWIFTGalaxy, and are therefore relatively memory-intensive. Masking individual cosmo_array datasets with __getitem__() avoids this: only masked copies of the individual arrays are returned in this case.

Parameters:

mask_collection (swiftgalaxy.masks.MaskCollection) – Set of masks to be applied to each particle type. Particle types may be omitted by setting their mask to None, or simply omitting them from the swiftgalaxy.masks.MaskCollection.

open_file() ContextManager[File]

Return a context manager that can be used to read the file.

This will use the existing handle if it is open. If not, we assume that we’re reading local HDF5 files using h5py and create a temporary handle.

Returns:

A context manager which can be used to read the file.

Return type:

ContextManager

recentre(new_centre: cosmo_array) None[source]

Set a new centre for the particle spatial coordinates.

The provided (spatial) coordinate centre is set to zero by subtracting it from the particle coordinates. Note that this is the new centre in the current coordinate system (not e.g. the simulation box coordinates). If the coordinate offset to be applied is known, use translate() instead. All datasets specified in the transforms_like_coordinates argument to SWIFTGalaxy are transformed (by default coordinates for all present particle types).

Parameters:

new_centre (cosmo_array) – The new centre for the (spatial) coordinate system.

See also

translate()

recentre_velocity(new_centre: cosmo_array) None[source]

Recentre the velocity coordinates.

The provided velocity coordinate is set to zero by subtracting it from the particle velocities. Note that this is the new velocity centre in the current coordinate system (not e.g. the simulation box coordinates). If the velocity offset to be applied is known, use boost() instead. All datasets specified in the transforms_like_velocities argument to SWIFTGalaxy are transformed (by default velocities for all present particle types).

Parameters:

new_centre (cosmo_array) – The new centre for the velocity coordinate system.

See also

boost()

rotate(rotation: Rotation) None[source]

Apply a rotation to the particle spatial coordinates.

The provided rotation is applied to all particle coordinates. All datasets specified in the transforms_like_coordinates and transforms_like_velocities arguments to SWIFTGalaxy are transformed (by default coordinates and velocities for all present particle types).

Parameters:

rotation (scipy.spatial.transform.Rotation) – The rotation to be applied. Rotation supports several input formats, including axis-angle, rotation matrices, and others.

property rotation: Rotation

The current rotation of the coordinate frame.

Returns:

The current rotation.

Return type:

scipy.spatial.transform.Rotation

snapshot_filename: str
transforms_like_coordinates: Set[str]
transforms_like_velocities: Set[str]
translate(translation: cosmo_array) None[source]

Apply a translation to the particle spatial coordinates.

The provided translation vector is added to all particle coordinates. If the new centre position that should be set to zero is known, use recentre() instead. All datasets specified in the transforms_like_coordinates argument to SWIFTGalaxy are transformed (by default coordinates for all present particle types).

Parameters:

translation (cosmo_array) – The vector to translate by.

See also

recentre()

velocities_dataset_name: str
property velocity_centre: cosmo_array

Get the current origin of the velocity reference frame.

It is given with respect to the native simulation velocity reference frame.

Returns:

The origin of the velocity reference frame.

Return type:

cosmo_array

wrap_box() None[source]

Wrap the particle coordinates in a periodic box.

Recentres a particle coordinates from a periodic simulation volume such that the coordinate (0, 0, 0) is in the centre and the axes of the volume are aligned with the coordinate axes.

Note

This is invoked automatically after any coordinate translations or rotations, so manually calling this function should usually not be necessary.

class swiftgalaxy.Standalone(centre: cosmo_array | None = None, velocity_centre: cosmo_array | None = None, spatial_offsets: cosmo_array | None = None, extra_mask: str | MaskCollection | None = None)[source]

Bases: _HaloCatalogue

Use to initialize a SWIFTGalaxy without a halo catalogue.

Provides an interface to specify the minimum required information to instantiate a SWIFTGalaxy.

Parameters:
  • centre (cosmo_array, default: None) – A value for this parameter is required, and must have units of length. Specifies the geometric centre in simulation coordinates. Particle coordinates will be shifted such that this position is located at (0, 0, 0) and if the boundary is periodic it will be wrapped to place the origin at the centre.

  • velocity_centre (cosmo_array, default: None) – A value for this parameter is required, and must have units of speed. Specifies the reference velocity relative to the simulation frame. Particle velocities will be shifted such that a particle with the specified velocity in the simulation frame will have zero velocity in the SWIFTGalaxy frame.

  • spatial_offsets (cosmo_array, default: None) – Offsets along each axis to select a spatial region around the centre. May be used in conjunction with extra_mask, for example to select all simulation particles in an aperture around the object of interest (see ‘Masking’ section of documentation for a cookbook example). Provide a pair of offsets from the object’s centre along each axis to define the region, for example for a cube extending +/- 1 Mpc from the centre: cosmo_array([[-1.0, 1.0], [-1.0, 1.0], [-1.0, 1.0]], u.Mpc).

  • extra_mask (str or MaskCollection (optional), default: None) – Mask to apply to particles after spatial masking. A user-defined mask can be provided as an an object (such as a swiftgalaxy.masks.MaskCollection) that has attributes with names corresponding to present particle names (e.g. gas, dark_matter, etc.), each containing a mask.

Examples

Often the most pragmatic way to create a selection of particles using Standalone is to first select a spatial region guaranteed to contain the particles of interest and then create the final mask programatically using SWIFTGalaxy’s masking features. For example, suppose we know that there is a galaxy with its centre at (2, 2, 2) Mpc and we eventually want all particles in a spherical aperture 1 Mpc in radius around this point. We start with a cubic spatial mask enclosing this region:

from swiftgalaxy import SWIFTGalaxy, Standalone, MaskCollection
from swiftsimio import cosmo_array
import unyt as u

sg = SWIFTGalaxy(
    "my_snapshot.hdf5",
    Standalone(
        centre=cosmo_array([2.0, 2.0, 2.0], u.Mpc),
        velocity_centre=cosmo_array([0.0, 0.0, 0.0], u.km / u.s),
        spatial_offsets=cosmo_array(
            [[-1.0, 1.0], [-1.0, 1.0], [-1.0, 1.0]],
            u.Mpc,
        ),
        extra_mask=None,  # we'll define the exact set of particles later
    )
)

We next define the masks selecting particles in the desired spherical aperture, conveniently using SWIFTGalaxy’s spherical coordinates feature, and store them in a MaskCollection:

mask_collection = MaskCollection(
    gas=sg.gas.spherical_coordinates.r < 1 * u.Mpc,
    dark_matter=sg.dark_matter.spherical_coordinates.r < 1 * u.Mpc,
    stars=sg.stars.spherical_coordinates.r < 1 * u.Mpc,
    black_holes=sg.black_holes.spherical_coordinates.r < 1 * u.Mpc,
)

Finally, we apply the mask to the sg object:

sg = sg.mask_particles(mask_collection)

We’re now ready to proceed with analysis of the particles in the 1 Mpc spherical aperture using this sg object.

Note

mask_particles() applies the masks in-place. The mask could also be applied with the __getattr__() method (i.e. in square brackets), but this returns a copy of the SWIFTGalaxy object. If memory efficiency is a concern, prefer the mask_particles() approach.

property centre: cosmo_array

Obtain the centre specified at initialization.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The centre(s) of the object(s) of interest.

Return type:

cosmo_array

property count: int

Number of galaxies in the catalogue.

When in multi-galaxy mode and not masked to a single row this can be >1.

Returns:

The number of galaxies in the catalogue.

Return type:

int

property velocity_centre: cosmo_array

Obtain the velocity centre specified at initialization.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The velocity coordinate origin.

Return type:

cosmo_array

class swiftgalaxy.Velociraptor(velociraptor_filebase: str | None = None, velociraptor_files: dict | None = None, halo_index: int | Sequence[int] | None = None, extra_mask: str | MaskCollection = 'bound_only', centre_type: str = 'minpot', custom_spatial_offsets: cosmo_array | None = None)[source]

Bases: _HaloCatalogue

Interface to velociraptor halo catalogues for use with swiftgalaxy.

Takes a set of velociraptor output files and configuration options and provides an interface that swiftgalaxy understands. Also exposes the halo/galaxy properties computed by velociraptor for a single object of interest with the API provided by the velociraptor python package. Reading of properties is done on-the-fly, and only rows corresponding to the object of interest are read from disk.

Parameters:
  • velociraptor_filebase (str) – The initial part of the velociraptor filenames (possibly including path), e.g. if there is a halos.properties file, pass halos as this argument. Provide this or velociraptor_files, not both.

  • velociraptor_files (dict) – A dictionary containing the names of the velociraptor files (possibly including paths). There should be two entries, with keys properties and catalog_groups containing locations of the {halos}.properties and {halos}.catalog_groups files, respectively. Provide this or velociraptor_filebase, not both.

  • halo_index (int or list) – Position(s) of the object(s) of interest in the catalogue arrays. In the case of multiple objects, duplicate entries are not allowed.

  • extra_mask (str or MaskCollection (optional), default: "bound_only") – Mask to apply to particles after spatial masking. If "bound_only", then the galaxy is masked to include only the gravitationally bound particles as determined by velociraptor. A user-defined mask can also be provided as an an object (such as a swiftgalaxy.masks.MaskCollection) that has attributes with names corresponding to present particle names (e.g. gas, dark_matter, etc.), each containing a mask.

  • centre_type (str (optional), default: "minpot") – Type of centre, chosen from those provided by velociraptor. Default is the position of the particle with the minimum potential, "minpot"; other possibilities may include "", "_gas", "_star", "mbp" (most bound particle).

  • custom_spatial_offsets (~swiftsimio.objects.cosmo_array (optional), default: None) – A region to override the automatically-determined region enclosing group member particles. May be used in conjunction with extra_mask, for example to select all simulation particles in an aperture around the object of interest (see ‘Masking’ section of documentation for a cookbook example). Provide a pair of offsets from the object’s centre along each axis to define the region, for example for a cube extending +/- 1 Mpc from the centre: cosmo_array([[-1.0, 1.0], [-1.0, 1.0], [-1.0, 1.0]], u.Mpc).

Notes

Note

velociraptor only supports index access to catalogue arrays, not identifier access. This means that the halo_index is simply the position of the object of interest in the catalogue arrays.

Examples

Given a file halos.properties (and also halos.catalog_groups, etc.) at /output/path/, the following creates a Velociraptor object for the entry at index 3 in the catalogue (i.e. the 4th row, indexed from 0) and demonstrates retrieving its virial mass.

>>> cat = Velociraptor(
>>>     velociraptor_filebase="/output/path/halos",  # halos.properties file is at
>>>                                                  # /output/path/
>>>     halo_index=3,  # 4th entry in catalogue (indexed from 0)
>>> )
>>> cat
Masked velociraptor catalogue at /path/to/output/out.properties.
Contains the following field collections: metallicity, ids, energies,
stellar_age, spherical_overdensities, rotational_support,
star_formation_rate, masses, eigenvectors, radii, temperature, veldisp,
structure_type, velocities, positions, concentration, rvmax_quantities,
angular_momentum, projected_apertures, apertures,
element_mass_fractions, dust_mass_fractions, number,
hydrogen_phase_fractions, black_hole_masses, stellar_birth_densities,
snii_thermal_feedback_densities, species_fractions,
gas_hydrogen_species_masses, gas_H_and_He_masses,
gas_diffuse_element_masses, dust_masses_from_table, dust_masses,
stellar_luminosities, cold_dense_gas_properties,
log_element_ratios_times_masses, lin_element_ratios_times_masses,
element_masses_in_stars, fail_all
>>> cat.masses
Contains the following fields: mass_200crit, mass_200crit_excl,
mass_200crit_excl_gas, mass_200crit_excl_gas_nsf,
mass_200crit_excl_gas_sf, mass_200crit_excl_star, mass_200crit_gas,
mass_200crit_gas_nsf, mass_200crit_gas_sf, mass_200crit_star,
mass_200mean, mass_200mean_excl, mass_200mean_excl_gas,
mass_200mean_excl_gas_nsf, mass_200mean_excl_gas_sf,
mass_200mean_excl_star, mass_200mean_gas, mass_200mean_gas_nsf,
mass_200mean_gas_sf, mass_200mean_star, mass_bn98, mass_bn98_excl,
mass_bn98_excl_gas, mass_bn98_excl_gas_nsf, mass_bn98_excl_gas_sf,
mass_bn98_excl_star, mass_bn98_gas, mass_bn98_gas_nsf,
mass_bn98_gas_sf, mass_bn98_star, mass_fof, mass_bh, mass_gas,
mass_gas_30kpc, mass_gas_500c, mass_gas_rvmax, mass_gas_hight_excl,
mass_gas_hight_incl, mass_gas_incl, mass_gas_nsf, mass_gas_nsf_incl,
mass_gas_sf, mass_gas_sf_incl, mass_star, mass_star_30kpc,
mass_star_500c, mass_star_rvmax, mass_star_incl, mass_tot,
mass_tot_incl, mvir
>>> cat.masses.mvir
unyt_array(14.73875777, '10000000000.0*Msun')
property centre: cosmo_array

Obtain the centre specified by the centre_type from the halo catalogue.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The centre(s) of the object(s) of interest.

Return type:

cosmo_array

centre_type: str
property count: int

Number of galaxies in the catalogue.

When in multi-galaxy mode and not masked to a single row this can be >1.

Returns:

The number of galaxies in the catalogue.

Return type:

int

property halo_index: List[int] | int

Get the index of the objects of interest in the halo catalogue.

In multi-galaxy mode when no mask is active this gets the list of indices.

This is just the position in the arrays stored in the Velociraptor catalogue files.

Returns:

The index or indices of the object(s) of interest in the halo catalogue.

Return type:

int or list

velociraptor_files: Dict[str, str]
property velocity_centre: cosmo_array

Obtain the velocity centre from the halo catalogue.

Specified by the velocity_centre_type.

In multi-galaxy mode if no mask is active return the centres of all objects of interest, otherwise return the centre of the (current) object of interest.

Returns:

The centre(s) of the object(s) of interest.

Return type:

cosmo_array

velocity_centre_type: str

Submodules