SimSpin.jl Documentation

Usage

For the installation procedure of the SimSpin package please follow the installation instructions on the package's README.

Mock observables of an N-body simulation

Once installed, only five steps are required to take an observation, generate a mock-IFU datacube and then export it to a FITS file:

  1. Construct a Telescope object. This specifies the field of view to be used, the aperture shape, etc. In this example we will use the default SAMI telescope constructor. See Telescope Constructors for other default and customisable constructors.

    > telescope = SAMI()

  2. Construct an environment in which the observation is taken. This specifies the redshift of the galaxy, its inclination, the virial radius, the mass to light ratio and the seeing conditions respectively. We will set redshift to be 0.05, inclination to be 70 degrees, virial radius to be 200 kpc, to mass to light ratio to be 1 and use no blurring. See Environment Constructor and Blur Constructors for more details.

    > environment = Environment(0.05, 70, 200, 1.)

  3. Read in a simulation's particle data. In this example we will use the example file embedded in SimSpin. Usually you will need to pass sim_data the file path to your HDF5 data file. See Data Import for more details.

    > data = sim_data()

  4. Build the datacube as a combination of the galaxy particle data, a telescope and an environment. It returns a mock-IFU datacube and a summary of the observational properties used. See build_datacube for more details.

    > datacube, observe = build_datacube(data, telescope, environment)

  5. Export the datacube to a FITS file for viewing. See Data Export for more details

    > sim_FITS(datacube, observe, "SimSpin_Example_Observation.fits")

Note

The most efficient way to take multiple observations is to pass build_datacube an array of Environments. These can be created individually or by passing the Environment Constructor an array of the values that you wish to observe. For example, to take an observation at redshifts 0.05, 0.1 and 0.15 each at an inclination of 70 degrees and 80 degrees using a constant virial radius of 200 kpc there are two options:

    > environments = Environment([0.05; 0.1; 015], [70; 80], 200)
    > environments = Environment(0.05:0.05:0.15, 70:10:80, 200)

Both options will return an array of environments which can then be used as the environment parameter in Step 4 above. This will then return an array of tuples, each one consisting of the datacube and the observational properties used. See build_datacube and Environment Constructor for more details.

Mock observables of a hydro-dynamical simulation

To use SimSpin with a hydro-dynamical simulation a similar process is followed with notable, minor differences:

  1. Construct a Telescope object. We must specify a filter when using a hydro-dynamical simulation. See Telescope Constructors for more options.

    > telescope = SAMI(filter="g")

  2. Construct an environment in which the observation is taken. This is the same as for an n-body simulation but does not require a mass to light ratio. See Environment Constructor and Blur Constructors for more details.

    > environment = Environment(0.05, 70, 200)

  3. Read in a simulation's particle data. For hydro-dynamical simulations we must specify that we want SSP to be used. The default is to not use SSP data so if SSP is not actively set to true a hydro-dynamical simulation will still be observed but a mass to light ratio will be used instead of spectra. See Data Import for more details.

    > data = sim_data("your/filename.hdf5", ssp=true)

  4. Build the datacube as a combination of the galaxy particle data, a telescope and an environment. It returns a mock-IFU datacube and a summary of the observational properties used. See build_datacube for more details.

    > datacube, observe = build_datacube(data, telescope, environment)

  5. Export the datacube to a FITS file for viewing. See Data Export for more details

    > sim_FITS(datacube, observe, "SimSpin_Example_Observation.fits")

Warning

Do not multithread any SimSpin functions as a user. Each observation is already multithreaded and will run on as many threads as available. To see how to allow SimSpin to use more threads see Multi-Threading.

General Functions

SimSpin.build_datacubeFunction
build_datacube(galaxy_data, ifu, envir)

Returns a simulated ifu datacube for input, galaxy_data, an array of Particle and a summary of observation properties used. The instrument to be used for observation is given in the Telescope class parameter and the environmental variables (redshift, inclination, seeing, etc) are given by the Environment class parameter.

Parameters:

galaxy_data         Array of `Particle` describing galaxy.
ifu                 Struct of type `Telescope`.
envir               Struct of type `Environment`.

Returns:

datacube           3D array of a simulated IFU datacube with spatial and velocity binned fluxes.
observe             A struct of type `Observation` summarising the observational properties used.
source
build_datacube(galaxy_data, ifu, environment_array)

Calls build_datacube using each environment in the array of Environment types. This is the most efficient way to take multiple observations.

Parameters:

galaxy_data         Array of `Particle` describing galaxy.
ifu                 Struct of type `Telescope`.
envir               Array of `Environment` describing each observation to be taken.

Returns:

datacubes   An array of tuples. Each tuple constists of a datacube [1] and an `Observation` [2].
source
SimSpin.flux_gridFunction
flux_grid(parts_in_cell,
            observe,
            filter)

Computes the fluxes for each element of the IFU data-cube.

The purpose of this function is to construct the mock flux values within each cell of the IFU cube. It accepts output parameters from obs_data_prep() and returns a 3D array containing the flux at each cell position due to contributions from the particles. If ssd particles are supplied, an SED is generated in each cell using ProSpect. Else, the luminosity in each cell is converted to flux.

Parameters:

parts_in_cell       1D array of the particles corresponding to each element in the IFU data-cube.
observe             Struct of type `Observation` containing all observation parameters.
filter_value        If ssp particles are supplied, the filter within which the SED is generated.
source
SimSpin.ifu_cubeFunction
ifu_cube(flux_grid,
            parts_in_cell,
            observe)

The purpose of this function is to construct an IFU data cube. It accepts a flux grid in the format output by the flux_grid() function and returns a similar, IFU-like, 3D array where each particle's flux contributes a Gaussian distribution in the velocity axis.

Parameters:

flux_grid       Flux grid output by `flux_grid()`
parts_in_cell   1D array of the particles corresponding to each element in the IFU data-cube.
observe         Struct of type `Observation` containing all observation parameters.
source
SimSpin.obs_data_prepFunction
obs_data_prep(galaxy_data, ifu, envir)

This function prepares the particle data for a given observation with the given telescope.

Parameters:

galaxy_data         Array of `Particle` describing galaxy
ifu                 Struct of type `Telescope`
envir               Struct of type `Environment`

Returns:

galaxy_data         Array of particle data formatted for the observation specified by the parameters.
parts_in_cell       3D array of the particles corresponding to each element in the IFU data-cube.
observe             Struct of type `Observation` containing all observation parameters.
source
SimSpin.sim_to_galaxyFunction
sim_to_galaxy(sim_data; centre)

Returns particle data wrt centre of galaxy in both spatial and velocity space instead of wrt arbitrary point in simulation. Centre in position or centre in both position and velocity can be optionally specified. Else, the median of each distribution will be used as centre.

source

Constructors

Telescope Constructors

An IFU object denotes all the parameters required to make a mock observation. Any generic IFU can be made using the IFU() function. Default constructors can also be used to emulate famous IFU survey instruments. Currently SAMI, MaNGA, CALIFA and Hector are supported.

SimSpin.IFUType
IFU(fov, ap_shape, central_wvl, lsf_fwhm, pixel_sscale, pixel_vscale,  filter)

Creates a customisable, mock IFU telescope which can be used to "observe" simulations.

Parameters:

fov             The field of view of the IFU, diameter in arcseconds.
ap_shape        The shape of the field of view, with options "circular", "square" or "hexagonal".
central_wvl     The central filter wavelength used for the observation, given in angstroms.
lsf_fwhm        The line spread function full-width half-max, given in angstroms.
pixel_sscale    The corresponding spatial pixel scale associated with a given telescope output in arcseconds.
pixel_vscale    The corresponding velocity pixel scale associated with a given telescope filter output in angstroms.
filter          Optional. If particles type is ssp, the filter within which the SED is generated. Options include "r" and "g"  for SDSS-r and SDSS-g bands respectively.
source
SimSpin.SAMIFunction
SAMI(;filter)

Creates an IFU using parameters of the SAMI survey. Optional filters "r" or "g" can also be specified for use with SSP particles.

source
SimSpin.MaNGAFunction
MaNGA(;filter)

Creates an IFU using parameters of the MaNGA survey. Optional filters "r" or "g" can also be specified for use with SSP particles.

source
SimSpin.CALIFAFunction
CALIFA(;filter)

Creates an IFU using parameters of the CALIFA survey. Optional filters "r" or "g" can also be specified for use with SSP particles.

source
SimSpin.HectorFunction
Hector(;filter)

Creates an IFU using parameters of the CALIFA survey. Optional filters "r" or "g" can also be specified for use with SSP particles.

source

Environment Constructor

SimSpin.EnvironmentType
Environment(z, inc_deg, r200, mass2light, blur)

Creates a struct containing environmental parameters required for a mock observation of a simulated galaxy. Each argument can either be an array or a single value. If an array is passed then a corresponding array of Environments will be returned containing every permutation of the given parameters.

Parameters:

z               The projected redshift at which the observation is made.
inc_deg         The inclination at which to observe the galaxy in degrees. Relative to face on, rotated around semi-major axis.
r200            The virial radius specified in the simulation, kpc.
mass2light      Optional. The mass to light ratio for non-ssp, luminous particles. Can be omitted, a single value or a tuple to specifiy disk and bulge values separately.
blur            Optional. Struct of type `Blur` containing seeing information. If ommitted no blurring is used.
source

Blur Constructors

SimSpin.Gaussian_blurType
Gaussian_blur(;sigma, fwhm)

Create a struct containing seeing information.

Keyword arguments (only one may be specified):

sigma       The standard deviation of the point spread function
fwhm        The full width half max of the point spread function
source
SimSpin.Moffat_blurType
Moffat_blur(β;
            α,
            fwhm)

Create a struct containing seeing information. β and either α or fwhm must be specified. If both α and fwhm are specified, α is prioritised.

Arguments:

β           The power component in the Moffat distribution
α           The core width of the Moffat distribution (optional)
fwhm        The full width half max of the Moffat distribution (optional)
source

Multi-Threading

The SimSpin package has multi-threading enabled in some critical functions. To use SimSpin with more than one thread the only thing you need to do is start Julia with the desired number of threads. Instructions on how to do this can be found here.

Warning

Do not manually multithread any SimSpin functions as a user. Each observation is already multithreaded and will run on as many threads as available.

Data Import

To read in a SimSpin HDF5 file as defined below the sim_data function must be used.

SimSpin.sim_dataFunction
sim_data(filename;
        pytpe = [],
        ssp = false)

Reads in a SimSpin format HDF5 file at location, filename. Returns array of Sim_particles.

Keyword arguments (optional):

ptype       A vector of the particles types to be read in e.g. ptype = [1,3].
            If omitted all particles types will be read.
ssp         Boolean value to use ssp particle information.
source
sim_data(;ptype::Vector{} = [],
            ssp::Bool = false)

if name filename is specified sim_data returns particle data from SimSpin's example file.

Keyword arguments (optional):

ptype       A vector of the particles types to be read in e.g. ptype = [1,3].
            If omitted all particles types will be read.
ssp         Boolean value to use ssp particle information.
source

The expected file format accepted by SimSpin is outlined below. If you would like to generate this file automatically, a short Python function has been written that uses the pynbody package to read in various simulation data types and generate a SimSpin compatible HDF5 file. See create_SimSpinFile.

If you would rather generate the SimSpin file independently, the expected file format is outlined below.

> SimSpin_example.hdf5

>> /PartType0           # Each particle type included in the simulation has its own group.
>>> /PartType0/Mass     # Each group then has a series of data sets assocaited,
>>> /PartType0/vx       #   including the position, velocity and Mass of each particle.
>>> /PartType0/vy
>>> /PartType0/vz
>>> /PartType0/x
>>> /PartType0/y
>>> /PartType0/z

>> /PartType1
>>> ...

We use the same PartType definition as Gadget: PartTypeX where 0 - gas, 1 - dark matter, 2 - disc, 3 - bulge, 4 - stars. For PartType0-3, each PartType group contains the same data sets as above. If the simulation contains stars, the Age and Metallicity information for each particle is also included:

> SimSpin_example.hdf5
>> /PartType4
>>> /PartType4/Age
>>> /PartType4/Mass
>>> /PartType4/Metallicity
>>> /PartType4/vx        
>>> /PartType4/vy
>>> /PartType4/vz
>>> /PartType4/x
>>> /PartType4/y
>>> /PartType4/z

Optionally SSP Star particles can also be defined with initial mass and stellar formation time:

> SimSpin_example.hdf5
>> /PartType4
>>> /PartType4/InitialMass
>>> /PartType4/Mass
>>> /PartType4/Metallicity
>>> /PartType4/StellarFormationTime
>>> /PartType4/vx        
>>> /PartType4/vy
>>> /PartType4/vz
>>> /PartType4/x
>>> /PartType4/y
>>> /PartType4/z

If the file is set up in this way, the simulation data can easily be read into the SimSpin package.

References

A. Pontzen, R Roskar, G. Stinson and R. Woods, (2013), "pynbody: N-Body/SPH analysis for python", Astrophysics Source Code Library, record ascl:1305.002

Data Export

Currently only FITS file exports are supported.

FITS Export

SimSpin.sim_FITSFunction
sim_FITS(data_cube,
            observe,
            out_file)

Outputs FITS file containing 3D, velocity binned datacube.

Parameters:

data_cube   3D array of a simulated IFU datacube with spatial and velocity binned fluxes
observe     A struct of type `Observation` summarising the observational properties used.
out_file    String denoting path and name of FITS file to be output.
source
sim_FITS(data_cube,
            out_file)

Outputs FITS file containing 3D, velocity binned datacube.

Parameters:

data_cube   Tuple with simulated IFU datacube and Observation as output from `build_datacube()`.
out_file    String denoting path and name of FITS file to be output.
source
sim_FITS(out_data;
            folder = "~",
            file_prefix = "SimSpin")

Can be used to export multiple datacubes in different FITS files.

Parameters:

out_data    Array of (datacube, observation) tuples. As output by `build_datacube()` when taking multiple observations.
folder      Optional. The file path to where the output should be. Defaults to "~".
file_prefix Optional. Prefix for the start of the .fits file's name.
source