trident.make_simple_ray

trident.make_simple_ray(dataset_file, start_position, end_position, lines=None, ftype='gas', fields=None, solution_filename=None, data_filename=None, trajectory=None, redshift=None, setup_function=None, load_kwargs=None, line_database=None, ionization_table=None)[source]

Create a yt LightRay object for a single dataset (eg CGM). This is a wrapper function around yt’s LightRay interface to reduce some of the complexity there.

A simple ray is a straight line passing through a single dataset where each gas cell intersected by the line is sampled for the desired fields and stored. Several additional fields are created and stored including dl to represent the path length in space for each element in the ray, v_los to represent the line of sight velocity along the ray, and redshift, redshift_dopp, and redshift_eff to represent the cosmological redshift, doppler redshift and effective redshift (combined doppler and cosmological) for each element of the ray.

A simple ray is typically specified by its start and end positions in the dataset volume. Because a simple ray only probes a single output, it lacks foreground absorbers between the observer at z=0 and the redshift of the dataset that one would naturally encounter. Thus it is usually only appropriate for studying the circumgalactic medium rather than the intergalactic medium.

This function can accept a yt dataset already loaded in memory, or it can load a dataset if you pass it the dataset’s filename and optionally any load_kwargs or setup_function necessary to load/process it properly before generating the LightRay object.

The :lines: keyword can be set to automatically add all fields to the resulting ray necessary for later use with the SpectrumGenerator class. If the necessary fields do not exist for your line of choice, they will be added to your dataset before adding them to the ray.

If using the :lines: keyword with an SPH dataset, it is very important to set the :ftype: keyword appropriately, or you may end up calculating ion fields by interpolating on data already smoothed to the grid. This is generally not desired.

Parameters

Dataset_file:

string or yt Dataset object

Either a yt dataset or the filename of a dataset on disk. If you are passing it a filename, consider usage of the load_kwargs and setup_function kwargs.

Start_position, end_position:
 

list of floats or YTArray object

The coordinates of the starting and ending position of the desired ray. If providing a raw list, coordinates are assumed to be in code length units, but if providing a YTArray, any units can be specified.

lines=None, fields=None, solution_filename=None, data_filename=None, trajectory=None, redshift=None, line_database=None, ftype=”gas”, setup_function=None, load_kwargs=None, ionization_table=None):

Lines:

list of strings, optional

List of strings that determine which fields will be added to the ray to support line deposition to an absorption line spectrum. List can include things like “C”, “O VI”, or “Mg II ####”, where #### would be the integer wavelength value of the desired line. If set to ‘all’, includes all possible ions from H to Zn. :lines: can be used in conjunction with :fields: as they will not override each other. If using the :lines: keyword with an SPH dataset, it is very important to set the :ftype: keyword appropriately, or you may end up calculating ion fields by interpolating on data already smoothed to the grid. This is generally not desired. Default: None

Ftype:

string, optional

For use with the :lines: keyword. It is the field type of the fields to be added. It is the first string in the field tuple e.g. “gas” in (“gas”, “O_p5_number_density”). For SPH datasets, it is important to set this to the field type of the gas particles in your dataset (e.g. ‘PartType0’), as it determines the source data for the ion fields to be added. If you leave it set to “gas”, it will calculate the ion fields based on the hydro fields already smoothed on the grid, which is usually not desired. Default: “gas”

Fields:

list of strings, optional

The list of which fields to store in the output LightRay. See :lines: keyword for additional functionality that will add fields necessary for creating absorption line spectra for certain line features. Default: None

Solution_filename:
 

string, optional

Output filename of text file containing trajectory of LightRay through the dataset. Default: None

Data_filename:

string, optional

Output filename for ray data stored as an HDF5 file. Note that at present, you must save a ray to disk in order for it to be returned by this function. If set to None, defaults to ‘ray.h5’. Default: None

Trajectory:

list of floats, optional

The (r, theta, phi) direction of the LightRay. Use either end_position or trajectory, but not both. Default: None

Redshift:

float, optional

Sets the highest cosmological redshift of the ray. By default, it will use the cosmological redshift of the dataset, if set, and if not set, it will use a redshift of 0. Default: None

Setup_function:

function, optional

A function that will be called on the dataset as it is loaded but before the LightRay is generated. Very useful for adding derived fields and other manipulations of the dataset prior to LightRay creation. Default: None

Load_kwargs:

dict, optional

Dictionary of kwargs to be passed to the yt “load” function prior to creating the LightRay. Very useful for many frontends like Gadget, Tipsy, etc. for passing in “bounding_box”, “unit_base”, etc. Default: None

Line_database:

string, optional

For use with the :lines: keyword. If you want to limit the available ion fields to be added to those available in a particular subset, you can use a LineDatabase. This means when you set :lines:=’all’, it will only use those ions present in the corresponding LineDatabase. If :LineDatabase: is set to None, and :lines:=’all’, it will add every ion of every element up to Zinc. Default: None

Ionization_table:
 

string, optional

For use with the :lines: keyword. Path to an appropriately formatted HDF5 table that can be used to compute the ion fraction as a function of density, temperature, metallicity, and redshift. When set to None, it uses the table specified in ~/.trident/config Default: None

Example

Generate a simple ray passing from the lower left corner to the upper right corner through some Gizmo dataset where gas particles are ftype=’PartType0’:

>>> import trident
>>> import yt
>>> ds = yt.load('path/to/dataset')
>>> ray = trident.make_simple_ray(ds,
... start_position=ds.domain_left_edge, end_position=ds.domain_right_edge,
... lines=['H', 'O', 'Mg II'], ftype='PartType0')