yt.frontends.enzo.simulation_handling module

class yt.frontends.enzo.simulation_handling.EnzoCosmology(hubble_constant, omega_matter, omega_lambda, omega_radiation, omega_curvature, initial_redshift, unit_registry=None)[source]

Bases: Cosmology

a_from_t(t)

Compute the scale factor for a given age of the Universe.

Parameters:

t (YTQuantity or float) – Time since the Big Bang. If a float is given, units are assumed to be seconds.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.a_from_t(4.0e17))
age_integrand(z)
angular_diameter_distance(z_i, z_f)

Following Hogg (1999), the angular diameter distance is ‘the ratio of an object’s physical transverse size to its angular size in radians.’

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.angular_diameter_distance(0.0, 1.0).in_units("Mpc"))
angular_scale(z_i, z_f)

The proper transverse distance between two points at redshift z_f observed at redshift z_i per unit of angular separation.

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.angular_scale(0.0, 1.0).in_units("kpc / arcsec"))
property arr
comoving_radial_distance(z_i, z_f)

The comoving distance along the line of sight to on object at redshift, z_f, viewed at a redshift, z_i.

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.comoving_radial_distance(0.0, 1.0).in_units("Mpccm"))
comoving_transverse_distance(z_i, z_f)

When multiplied by some angle, the distance between two objects observed at redshift, z_f, with an angular separation given by that angle, viewed by an observer at redshift, z_i (Hogg 1999).

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.comoving_transverse_distance(0.0, 1.0).in_units("Mpccm"))
comoving_volume(z_i, z_f)

“The comoving volume is the volume measure in which number densities of non-evolving objects locked into Hubble flow are constant with redshift.” – Hogg (1999)

Parameters:
  • z_i (float) – The lower redshift of the interval.

  • z_f (float) – The higher redshift of the interval.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.comoving_volume(0.0, 1.0).in_units("Gpccm**3"))
critical_density(z)

The density required for closure of the Universe at a given redshift in the proper frame.

Parameters:

z (float) – Redshift.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.critical_density(0.0).in_units("g/cm**3"))
>>> print(co.critical_density(0).in_units("Msun/Mpc**3"))
expansion_factor(z)

The ratio between the Hubble parameter at a given redshift and redshift zero.

This is also the primary function integrated to calculate the cosmological distances.

get_dark_factor(z)

This function computes the additional term that enters the expansion factor when using non-standard dark energy. See Dolag et al 2004 eq. 7 for ref (but note that there’s a typo in his eq. There should be no negative sign).

At the moment, this only works using the parameterization given in Linder 2002 eq. 7: w(a) = w0 + wa(1 - a) = w0 + wa * z / (1+z). This gives rise to an analytic expression. It is also only functional for Gadget simulations, at the moment.

Parameters:

z (float) – Redshift

hubble_distance()

The distance corresponding to c / h, where c is the speed of light and h is the Hubble parameter in units of 1 / time.

hubble_parameter(z)

The value of the Hubble parameter at a given redshift.

Parameters:

z (float) – Redshift.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.hubble_parameter(1.0).in_units("km/s/Mpc"))
inverse_expansion_factor(z)
lookback_time(z_i, z_f)

The difference in the age of the Universe between the redshift interval z_i to z_f.

Parameters:
  • z_i (float) – The lower redshift of the interval.

  • z_f (float) – The higher redshift of the interval.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.lookback_time(0.0, 1.0).in_units("Gyr"))
luminosity_distance(z_i, z_f)

The distance that would be inferred from the inverse-square law of light and the measured flux and luminosity of the observed object.

Parameters:
  • z_i (float) – The redshift of the observer.

  • z_f (float) – The redshift of the observed object.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.luminosity_distance(0.0, 1.0).in_units("Mpc"))
path_length(z_i, z_f)
path_length_function(z)
property quan
t_from_a(a)

Compute the age of the Universe for a given scale factor.

Parameters:

a (float) – Scale factor.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.t_from_a(1.0).in_units("Gyr"))
t_from_z(z)

Compute the age of the Universe for a given redshift.

Parameters:

z (float) – Redshift.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.t_from_z(0.0).in_units("Gyr"))
z_from_t(t)

Compute the redshift for a given age of the Universe.

Parameters:

t (YTQuantity or float) – Time since the Big Bang. If a float is given, units are assumed to be seconds.

Examples

>>> from yt.utilities.cosmology import Cosmology
>>> co = Cosmology()
>>> print(co.z_from_t(4.0e17))
class yt.frontends.enzo.simulation_handling.EnzoSimulation(outputs, *args, **kwargs)[source]

Bases: SimulationTimeSeries

Initialize an Enzo Simulation object.

Upon creation, the parameter file is parsed and the time and redshift are calculated and stored in all_outputs. A time units dictionary is instantiated to allow for time outputs to be requested with physical time units. The get_time_series can be used to generate a DatasetSeries object.

parameter_filenamestr

The simulation parameter file.

find_outputsbool

If True, subdirectories within the GlobalDir directory are searched one by one for datasets. Time and redshift information are gathered by temporarily instantiating each dataset. This can be used when simulation data was created in a non-standard way, making it difficult to guess the corresponding time and redshift information. Default: False.

Examples

>>> import yt
>>> es = yt.load_simulation("enzo_tiny_cosmology/32Mpc_32.enzo", "Enzo")
>>> es.get_time_series()
>>> for ds in es:
...     print(ds.current_time)
property arr
eval(tasks, obj=None)
classmethod from_output_log(output_log, line_prefix='DATASET WRITTEN', parallel=True)
get_time_series(time_data=True, redshift_data=True, initial_time=None, final_time=None, initial_redshift=None, final_redshift=None, initial_cycle=None, final_cycle=None, times=None, redshifts=None, tolerance=None, parallel=True, setup_function=None)[source]

Instantiate a DatasetSeries object for a set of outputs.

If no additional keywords given, a DatasetSeries object will be created with all potential datasets created by the simulation.

Outputs can be gather by specifying a time or redshift range (or combination of time and redshift), with a specific list of times or redshifts, a range of cycle numbers (for cycle based output), or by simply searching all subdirectories within the simulation directory.

time_databool

Whether or not to include time outputs when gathering datasets for time series. Default: True.

redshift_databool

Whether or not to include redshift outputs when gathering datasets for time series. Default: True.

initial_timetuple of type (float, str)

The earliest time for outputs to be included. This should be given as the value and the string representation of the units. For example, (5.0, “Gyr”). If None, the initial time of the simulation is used. This can be used in combination with either final_time or final_redshift. Default: None.

final_timetuple of type (float, str)

The latest time for outputs to be included. This should be given as the value and the string representation of the units. For example, (13.7, “Gyr”). If None, the final time of the simulation is used. This can be used in combination with either initial_time or initial_redshift. Default: None.

timestuple of type (float array, str)

A list of times for which outputs will be found and the units of those values. For example, ([0, 1, 2, 3], “s”). Default: None.

initial_redshiftfloat

The earliest redshift for outputs to be included. If None, the initial redshift of the simulation is used. This can be used in combination with either final_time or final_redshift. Default: None.

final_redshiftfloat

The latest redshift for outputs to be included. If None, the final redshift of the simulation is used. This can be used in combination with either initial_time or initial_redshift. Default: None.

redshiftsarray_like

A list of redshifts for which outputs will be found. Default: None.

initial_cyclefloat

The earliest cycle for outputs to be included. If None, the initial cycle of the simulation is used. This can only be used with final_cycle. Default: None.

final_cyclefloat

The latest cycle for outputs to be included. If None, the final cycle of the simulation is used. This can only be used in combination with initial_cycle. Default: None.

tolerancefloat

Used in combination with “times” or “redshifts” keywords, this is the tolerance within which outputs are accepted given the requested times or redshifts. If None, the nearest output is always taken. Default: None.

parallelbool/int

If True, the generated DatasetSeries will divide the work such that a single processor works on each dataset. If an integer is supplied, the work will be divided into that number of jobs. Default: True.

setup_functioncallable, accepts a ds

This function will be called whenever a dataset is loaded.

Examples

>>> import yt
>>> es = yt.load_simulation("enzo_tiny_cosmology/32Mpc_32.enzo", "Enzo")
>>> es.get_time_series(
...     initial_redshift=10, final_time=(13.7, "Gyr"), redshift_data=False
... )
>>> for ds in es:
...     print(ds.current_time)
>>> es.get_time_series(redshifts=[3, 2, 1, 0])
>>> for ds in es:
...     print(ds.current_time)
property outputs
particle_trajectories(indices, fields=None, suppress_logging=False, ptype=None)

Create a collection of particle trajectories in time over a series of datasets.

Parameters:
  • indices (array_like) – An integer array of particle indices whose trajectories we want to track. If they are not sorted they will be sorted.

  • fields (list of strings, optional) – A set of fields that is retrieved when the trajectory collection is instantiated. Default: None (will default to the fields ‘particle_position_x’, ‘particle_position_y’, ‘particle_position_z’)

  • suppress_logging (boolean) – Suppress yt’s logging when iterating over the simulation time series. Default: False

  • ptype (str, optional) – Only use this particle type. Default: None, which uses all particle type.

Examples

>>> my_fns = glob.glob("orbit_hdf5_chk_00[0-9][0-9]")
>>> my_fns.sort()
>>> fields = [
...     ("all", "particle_position_x"),
...     ("all", "particle_position_y"),
...     ("all", "particle_position_z"),
...     ("all", "particle_velocity_x"),
...     ("all", "particle_velocity_y"),
...     ("all", "particle_velocity_z"),
... ]
>>> ds = load(my_fns[0])
>>> init_sphere = ds.sphere(ds.domain_center, (0.5, "unitary"))
>>> indices = init_sphere[("all", "particle_index")].astype("int64")
>>> ts = DatasetSeries(my_fns)
>>> trajs = ts.particle_trajectories(indices, fields=fields)
>>> for t in trajs:
...     print(
...         t[("all", "particle_velocity_x")].max(),
...         t[("all", "particle_velocity_x")].min(),
...     )

Notes

This function will fail if there are duplicate particle ids or if some of the particle disappear.

piter(storage=None, dynamic=False)

Iterate over time series components in parallel.

This allows you to iterate over a time series while dispatching individual components of that time series to different processors or processor groups. If the parallelism strategy was set to be multi-processor (by “parallel = N” where N is an integer when the DatasetSeries was created) this will issue each dataset to an N-processor group. For instance, this would allow you to start a 1024 processor job, loading up 100 datasets in a time series and creating 8 processor groups of 128 processors each, each of which would be assigned a different dataset. This could be accomplished as shown in the examples below. The storage option is as seen in parallel_objects() which is a mechanism for storing results of analysis on an individual dataset and then combining the results at the end, so that the entire set of processors have access to those results.

Note that supplying a store changes the iteration mechanism; see below.

Parameters:
  • storage (dict) – This is a dictionary, which will be filled with results during the course of the iteration. The keys will be the dataset indices and the values will be whatever is assigned to the result attribute on the storage during iteration.

  • dynamic (boolean) – This governs whether or not dynamic load balancing will be enabled. This requires one dedicated processor; if this is enabled with a set of 128 processors available, only 127 will be available to iterate over objects as one will be load balancing the rest.

Examples

Here is an example of iteration when the results do not need to be stored. One processor will be assigned to each dataset.

>>> ts = DatasetSeries("DD*/DD*.index")
>>> for ds in ts.piter():
...     SlicePlot(ds, "x", ("gas", "density")).save()
...

This demonstrates how one might store results:

>>> def print_time(ds):
...     print(ds.current_time)
...
>>> ts = DatasetSeries("DD*/DD*.index", setup_function=print_time)
...
>>> my_storage = {}
>>> for sto, ds in ts.piter(storage=my_storage):
...     v, c = ds.find_max(("gas", "density"))
...     sto.result = (v, c)
...
>>> for i, (v, c) in sorted(my_storage.items()):
...     print("% 4i  %0.3e" % (i, v))
...

This shows how to dispatch 4 processors to each dataset:

>>> ts = DatasetSeries("DD*/DD*.index", parallel=4)
>>> for ds in ts.piter():
...     ProjectionPlot(ds, "x", ("gas", "density")).save()
...
print_key_parameters()

Print out some key parameters for the simulation.

property quan