yt.data_objects.selection_objects.region module¶
- class yt.data_objects.selection_objects.region.YTRegion(center, left_edge, right_edge, fields=None, ds=None, field_parameters=None, data_source=None)[source]¶
Bases:
YTSelectionContainer3D
A 3D region of data with an arbitrary center.
Takes an array of three left_edge coordinates, three right_edge coordinates, and a center that can be anywhere in the domain. If the selected region extends past the edges of the domain, no data will be found there, though the object’s left_edge or right_edge are not modified.
- Parameters:
center (array_like) – The center of the region
left_edge (array_like) – The left edge of the region
right_edge (array_like) – The right edge of the region
- apply_units(arr, units)¶
- argmax(field, axis=None)¶
Return the values at which the field is maximized.
This will, in a parallel-aware fashion, find the maximum value and then return to you the values at that maximum location that are requested for “axis”. By default it will return the spatial positions (in the natural coordinate system), but it can be any field
- Parameters:
field (string or tuple field name) – The field to maximize.
axis (string or list of strings, optional) – If supplied, the fields to sample along; if not supplied, defaults to the coordinate fields. This can be the name of the coordinate fields (i.e., ‘x’, ‘y’, ‘z’) or a list of fields, but cannot be 0, 1, 2.
- Return type:
A list of YTQuantities as specified by the axis argument.
Examples
>>> temp_at_max_rho = reg.argmax( ... ("gas", "density"), axis=("gas", "temperature") ... ) >>> max_rho_xyz = reg.argmax(("gas", "density")) >>> t_mrho, v_mrho = reg.argmax( ... ("gas", "density"), ... axis=[("gas", "temperature"), ("gas", "velocity_magnitude")], ... ) >>> x, y, z = reg.argmax(("gas", "density"))
- argmin(field, axis=None)¶
Return the values at which the field is minimized.
This will, in a parallel-aware fashion, find the minimum value and then return to you the values at that minimum location that are requested for “axis”. By default it will return the spatial positions (in the natural coordinate system), but it can be any field
- Parameters:
field (string or tuple field name) – The field to minimize.
axis (string or list of strings, optional) – If supplied, the fields to sample along; if not supplied, defaults to the coordinate fields. This can be the name of the coordinate fields (i.e., ‘x’, ‘y’, ‘z’) or a list of fields, but cannot be 0, 1, 2.
- Return type:
A list of YTQuantities as specified by the axis argument.
Examples
>>> temp_at_min_rho = reg.argmin( ... ("gas", "density"), axis=("gas", "temperature") ... ) >>> min_rho_xyz = reg.argmin(("gas", "density")) >>> t_mrho, v_mrho = reg.argmin( ... ("gas", "density"), ... axis=[("gas", "temperature"), ("gas", "velocity_magnitude")], ... ) >>> x, y, z = reg.argmin(("gas", "density"))
- property blocks¶
- calculate_isocontour_flux(field, value, field_x, field_y, field_z, fluxing_field=None)¶
This identifies isocontours on a cell-by-cell basis, with no consideration of global connectedness, and calculates the flux over those contours.
This function will conduct marching cubes on all the cells in a given data container (grid-by-grid), and then for each identified triangular segment of an isocontour in a given cell, calculate the gradient (i.e., normal) in the isocontoured field, interpolate the local value of the “fluxing” field, the area of the triangle, and then return:
area * local_flux_value * (n dot v)
Where area, local_value, and the vector v are interpolated at the barycenter (weighted by the vertex values) of the triangle. Note that this specifically allows for the field fluxing across the surface to be different from the field being contoured. If the fluxing_field is not specified, it is assumed to be 1.0 everywhere, and the raw flux with no local-weighting is returned.
Additionally, the returned flux is defined as flux into the surface, not flux out of the surface.
- Parameters:
field (string) – Any field that can be obtained in a data object. This is the field which will be isocontoured and used as the “local_value” in the flux equation.
value (float) – The value at which the isocontour should be calculated.
field_x (string) – The x-component field
field_y (string) – The y-component field
field_z (string) – The z-component field
fluxing_field (string, optional) – The field whose passage over the surface is of interest. If not specified, assumed to be 1.0 everywhere.
- Returns:
flux – The summed flux. Note that it is not currently scaled; this is simply the code-unit area times the fields.
- Return type:
Examples
This will create a data object, find a nice value in the center, and calculate the metal flux over it.
>>> dd = ds.all_data() >>> rho = dd.quantities["WeightedAverageQuantity"]( ... ("gas", "density"), weight=("gas", "cell_mass") ... ) >>> flux = dd.calculate_isocontour_flux( ... ("gas", "density"), ... rho, ... ("gas", "velocity_x"), ... ("gas", "velocity_y"), ... ("gas", "velocity_z"), ... ("gas", "metallicity"), ... )
- chunks(fields, chunking_style, **kwargs)¶
- clear_data()¶
Clears out all data from the YTDataContainer instance, freeing memory.
- clone()¶
Clone a data object.
This will make a duplicate of a data object; note that the field_parameters may not necessarily be deeply-copied. If you modify the field parameters in-place, it may or may not be shared between the objects, depending on the type of object that that particular field parameter is.
Notes
One use case for this is to have multiple identical data objects that are being chunked over in different orders.
Examples
>>> ds = yt.load("IsolatedGalaxy/galaxy0030/galaxy0030") >>> sp = ds.sphere("c", 0.1) >>> sp_clone = sp.clone() >>> sp["gas", "density"] >>> print(sp.field_data.keys()) [("gas", "density")] >>> print(sp_clone.field_data.keys()) []
- comm = None¶
- create_firefly_object(datadir=None, fields_to_include=None, fields_units=None, default_decimation_factor=100, velocity_units='km/s', coordinate_units='kpc', show_unused_fields=0, *, JSONdir=None, match_any_particle_types=True, **kwargs)¶
This function links a region of data stored in a yt dataset to the Python frontend API for [Firefly](http://github.com/ageller/Firefly), a browser-based particle visualization tool.
- Parameters:
datadir (string) – Path to where any .json files should be saved. If a relative path will assume relative to ${HOME}. A value of None will default to ${HOME}/Data.
fields_to_include (array_like of strings or field tuples) – A list of fields that you want to include in your Firefly visualization for on-the-fly filtering and colormapping.
default_decimation_factor (integer) – The factor by which you want to decimate each particle group by (e.g. if there are 1e7 total particles in your simulation you might want to set this to 100 at first). Randomly samples your data like shuffled_data[::decimation_factor] so as to not overtax a system. This is adjustable on a per particle group basis by changing the returned reader’s reader.particleGroup[i].decimation_factor before calling reader.writeToDisk().
velocity_units (string) – The units that the velocity should be converted to in order to show streamlines in Firefly. Defaults to km/s.
coordinate_units (string) – The units that the coordinates should be converted to. Defaults to kpc.
show_unused_fields (boolean) – A flag to optionally print the fields that are available, in the dataset but were not explicitly requested to be tracked.
match_any_particle_types (boolean) – If True, when any of the fields_to_include match multiple particle groups then the field will be added for all matching particle groups. If False, an error is raised when encountering an ambiguous field. Default is True.
to (Any additional keyword arguments are passed)
firefly.data_reader.Reader.__init__
- Returns:
reader – A reader object from the Firefly, configured to output the current region selected
- Return type:
Firefly.data_reader.Reader object
Examples
>>> ramses_ds = yt.load( ... "/Users/agurvich/Desktop/yt_workshop/" ... + "DICEGalaxyDisk_nonCosmological/output_00002/info_00002.txt" ... )
>>> region = ramses_ds.sphere(ramses_ds.domain_center, (1000, "kpc"))
>>> reader = region.create_firefly_object( ... "IsoGalaxyRamses", ... fields_to_include=[ ... "particle_extra_field_1", ... "particle_extra_field_2", ... ], ... fields_units=["dimensionless", "dimensionless"], ... )
>>> reader.settings["color"]["io"] = [1, 1, 0, 1] >>> reader.particleGroups[0].decimation_factor = 100 >>> reader.writeToDisk()
- cut_region(field_cuts, field_parameters=None, locals=None)¶
Return a YTCutRegion, where the a cell is identified as being inside the cut region based on the value of one or more fields. Note that in previous versions of yt the name ‘grid’ was used to represent the data object used to construct the field cut, as of yt 3.0, this has been changed to ‘obj’.
- Parameters:
field_cuts (list of strings) – A list of conditionals that will be evaluated. In the namespace available, these conditionals will have access to ‘obj’ which is a data object of unknown shape, and they must generate a boolean array. For instance, conditionals = [“obj[‘gas’, ‘temperature’] < 1e3”]
field_parameters (dictionary) – A dictionary of field parameters to be used when applying the field cuts.
locals (dictionary) – A dictionary of local variables to use when defining the cut region.
Examples
To find the total mass of hot gas with temperature greater than 10^6 K in your volume:
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.cut_region(["obj['gas', 'temperature'] > 1e6"]) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- exclude_above(field, value, units=None)¶
This function will return a YTCutRegion where all of the regions whose field is above a given value are masked.
- Parameters:
field (string) – The field in which the conditional will be applied.
value (float) – The minimum value that will not be masked in the output YTCutRegion.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field above the given value masked.
- Return type:
Examples
To find the total mass of hot gas with temperature colder than 10^6 K in your volume:
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_above(("gas", "temperature"), 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- exclude_below(field, value, units=None)¶
This function will return a YTCutRegion where all of the regions whose field is below a given value are masked.
- Parameters:
field (string) – The field in which the conditional will be applied.
value (float) – The minimum value that will not be masked in the output YTCutRegion.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field below the given value masked.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_below(("gas", "temperature"), 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- exclude_equal(field, value, units=None)¶
This function will return a YTCutRegion where all of the regions whose field are equal to given value are masked.
- Parameters:
field (string) – The field in which the conditional will be applied.
value (float) – The minimum value that will not be masked in the output YTCutRegion.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field equal to the given value masked.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_equal(("gas", "temperature"), 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- exclude_inside(field, min_value, max_value, units=None)¶
This function will return a YTCutRegion where all of the regions whose field are inside the interval from min_value to max_value.
- Parameters:
field (string) – The field in which the conditional will be applied.
min_value (float) – The minimum value inside the interval to be excluded.
max_value (float) – The maximum value inside the interval to be excluded.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field inside the given interval excluded.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_inside(("gas", "temperature"), 1e5, 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- exclude_nan(field, units=None)¶
This function will return a YTCutRegion where all of the regions whose field is NaN are masked.
- Parameters:
field (string) – The field in which the conditional will be applied.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the NaN entries of the field masked.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_nan(("gas", "temperature")) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- exclude_outside(field, min_value, max_value, units=None)¶
This function will return a YTCutRegion where all of the regions whose field are outside the interval from min_value to max_value.
- Parameters:
field (string) – The field in which the conditional will be applied.
min_value (float) – The minimum value inside the interval to be excluded.
max_value (float) – The maximum value inside the interval to be excluded.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field outside the given interval excluded.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_outside(("gas", "temperature"), 1e5, 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- extract_connected_sets(field, num_levels, min_val, max_val, log_space=True, cumulative=True)¶
This function will create a set of contour objects, defined by having connected cell structures, which can then be studied and used to ‘paint’ their source grids, thus enabling them to be plotted.
Note that this function can return a connected set object that has no member values.
- extract_isocontours(field, value, filename=None, rescale=False, sample_values=None)¶
This identifies isocontours on a cell-by-cell basis, with no consideration of global connectedness, and returns the vertices of the Triangles in that isocontour.
This function simply returns the vertices of all the triangles calculated by the marching cubes algorithm; for more complex operations, such as identifying connected sets of cells above a given threshold, see the extract_connected_sets function. This is more useful for calculating, for instance, total isocontour area, or visualizing in an external program (such as MeshLab.)
- Parameters:
field (string) – Any field that can be obtained in a data object. This is the field which will be isocontoured.
value (float) – The value at which the isocontour should be calculated.
filename (string, optional) – If supplied, this file will be filled with the vertices in .obj format. Suitable for loading into meshlab.
rescale (bool, optional) – If true, the vertices will be rescaled within their min/max.
sample_values (string, optional) – Any field whose value should be extracted at the center of each triangle.
- Returns:
verts (array of floats) – The array of vertices, x,y,z. Taken in threes, these are the triangle vertices.
samples (array of floats) – If sample_values is specified, this will be returned and will contain the values of the field specified at the center of each triangle.
Examples
This will create a data object, find a nice value in the center, and output the vertices to “triangles.obj” after rescaling them.
>>> dd = ds.all_data() >>> rho = dd.quantities["WeightedAverageQuantity"]( ... ("gas", "density"), weight=("gas", "cell_mass") ... ) >>> verts = dd.extract_isocontours( ... ("gas", "density"), rho, "triangles.obj", True ... )
- property fcoords¶
- property fcoords_vertex¶
- property fwidth¶
- get_data(fields=None)¶
- get_dependencies(fields)¶
- get_field_parameter(name, default=None)¶
This is typically only used by derived field functions, but it returns parameters used to generate fields.
- has_field_parameter(name)¶
Checks if a field parameter is set.
- has_key(key)¶
Checks if a data field already exists.
- property icoords¶
- include_above(field, value, units=None)¶
This function will return a YTCutRegion where only the regions whose field is above a given value are included.
- Parameters:
field (string) – The field in which the conditional will be applied.
value (float) – The minimum value that will not be masked in the output YTCutRegion.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field above the given value masked.
- Return type:
Examples
To find the total mass of hot gas with temperature warmer than 10^6 K in your volume:
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.include_above(("gas", "temperature"), 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- include_below(field, value, units=None)¶
This function will return a YTCutRegion where only the regions whose field is below a given value are included.
- Parameters:
field (string) – The field in which the conditional will be applied.
value (float) – The minimum value that will not be masked in the output YTCutRegion.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with only regions with the field below the given value included.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.include_below(("gas", "temperature"), 1e5, 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- include_equal(field, value, units=None)¶
This function will return a YTCutRegion where only the regions whose field are equal to given value are included.
- Parameters:
field (string) – The field in which the conditional will be applied.
value (float) – The minimum value that will not be masked in the output YTCutRegion.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field equal to the given value included.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.include_equal(("gas", "temperature"), 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- include_inside(field, min_value, max_value, units=None)¶
This function will return a YTCutRegion where only the regions whose field are inside the interval from min_value to max_value are included.
- Parameters:
field (string) – The field in which the conditional will be applied.
min_value (float) – The minimum value inside the interval to be excluded.
max_value (float) – The maximum value inside the interval to be excluded.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field inside the given interval excluded.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.include_inside(("gas", "temperature"), 1e5, 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- include_outside(field, min_value, max_value, units=None)¶
This function will return a YTCutRegion where only the regions whose field are outside the interval from min_value to max_value are included.
- Parameters:
field (string) – The field in which the conditional will be applied.
min_value (float) – The minimum value inside the interval to be excluded.
max_value (float) – The maximum value inside the interval to be excluded.
units (string or None) – The units of the value threshold. None will use the default units given in the field.
- Returns:
cut_region – The YTCutRegion with the field outside the given interval excluded.
- Return type:
Examples
>>> ds = yt.load("RedshiftOutput0005") >>> ad = ds.all_data() >>> cr = ad.exclude_outside(("gas", "temperature"), 1e5, 1e6) >>> print(cr.quantities.total_quantity(("gas", "cell_mass")).in_units("Msun"))
- property index¶
- integrate(field, weight=None, axis=None, *, moment=1)¶
Compute the integral (projection) of a field along an axis.
This projects a field along an axis.
- Parameters:
field (string or tuple field name) – The field to project.
weight (string or tuple field name) – The field to weight the projection by
axis (string) – The axis to project along.
moment (integer, optional) – for a weighted projection, moment = 1 (the default) corresponds to a weighted average. moment = 2 corresponds to a weighted standard deviation.
- Return type:
YTProjection
Examples
>>> column_density = reg.integrate(("gas", "density"), axis=("index", "z"))
- property ires¶
- keys()¶
- max(field, axis=None)¶
Compute the maximum of a field, optionally along an axis.
This will, in a parallel-aware fashion, compute the maximum of the given field. Supplying an axis will result in a return value of a YTProjection, with method ‘max’ for maximum intensity. If the max has already been requested, it will use the cached extrema value.
- Parameters:
field (string or tuple field name) – The field to maximize.
axis (string, optional) – If supplied, the axis to project the maximum along.
- Return type:
Either a scalar or a YTProjection.
Examples
>>> max_temp = reg.max(("gas", "temperature")) >>> max_temp_proj = reg.max(("gas", "temperature"), axis=("index", "x"))
- property max_level¶
- mean(field, axis=None, weight=None)¶
Compute the mean of a field, optionally along an axis, with a weight.
This will, in a parallel-aware fashion, compute the mean of the given field. If an axis is supplied, it will return a projection, where the weight is also supplied. By default the weight field will be “ones” or “particle_ones”, depending on the field being averaged, resulting in an unweighted average.
- Parameters:
field (string or tuple field name) – The field to average.
axis (string, optional) – If supplied, the axis to compute the mean along (i.e., to project along)
weight (string, optional) – The field to use as a weight.
- Return type:
Scalar or YTProjection.
Examples
>>> avg_rho = reg.mean(("gas", "density"), weight="cell_volume") >>> rho_weighted_T = reg.mean( ... ("gas", "temperature"), axis=("index", "y"), weight=("gas", "density") ... )
- min(field, axis=None)¶
Compute the minimum of a field.
This will, in a parallel-aware fashion, compute the minimum of the given field. Supplying an axis will result in a return value of a YTProjection, with method ‘min’ for minimum intensity. If the min has already been requested, it will use the cached extrema value.
- Parameters:
field (string or tuple field name) – The field to minimize.
axis (string, optional) – If supplied, the axis to compute the minimum along.
- Return type:
Either a scalar or a YTProjection.
Examples
>>> min_temp = reg.min(("gas", "temperature")) >>> min_temp_proj = reg.min(("gas", "temperature"), axis=("index", "x"))
- property min_level¶
- partition_index_2d(axis)¶
- partition_index_3d(ds, padding=0.0, rank_ratio=1)¶
- partition_index_3d_bisection_list()¶
Returns an array that is used to drive _partition_index_3d_bisection, below.
- partition_region_3d(left_edge, right_edge, padding=0.0, rank_ratio=1)¶
Given a region, it subdivides it into smaller regions for parallel analysis.
- property pf¶
- profile(bin_fields, fields, n_bins=64, extrema=None, logs=None, units=None, weight_field=('gas', 'mass'), accumulation=False, fractional=False, deposition='ngp', *, override_bins=None)¶
Create a 1, 2, or 3D profile object from this data_source.
The dimensionality of the profile object is chosen by the number of fields given in the bin_fields argument. This simply calls
yt.data_objects.profiles.create_profile()
.- Parameters:
bin_fields (list of strings) – List of the binning fields for profiling.
fields (list of strings) – The fields to be profiled.
n_bins (int or list of ints) – The number of bins in each dimension. If None, 64 bins for each bin are used for each bin field. Default: 64.
extrema (dict of min, max tuples) – Minimum and maximum values of the bin_fields for the profiles. The keys correspond to the field names. Defaults to the extrema of the bin_fields of the dataset. If a units dict is provided, extrema are understood to be in the units specified in the dictionary.
logs (dict of boolean values) – Whether or not to log the bin_fields for the profiles. The keys correspond to the field names. Defaults to the take_log attribute of the field.
units (dict of strings) – The units of the fields in the profiles, including the bin_fields.
weight_field (str or tuple field identifier) – The weight field for computing weighted average for the profile values. If None, the profile values are sums of the data in each bin.
accumulation (bool or list of bools) – If True, the profile values for a bin n are the cumulative sum of all the values from bin 0 to n. If -True, the sum is reversed so that the value for bin n is the cumulative sum from bin N (total bins) to n. If the profile is 2D or 3D, a list of values can be given to control the summation in each dimension independently. Default: False.
fractional (If True the profile values are divided by the sum of all) – the profile data such that the profile represents a probability distribution function.
deposition (Controls the type of deposition used for ParticlePhasePlots.) – Valid choices are ‘ngp’ and ‘cic’. Default is ‘ngp’. This parameter is ignored if the input fields are not of particle type.
override_bins (dict of bins to profile plot with) – If set, ignores n_bins and extrema settings and uses the supplied bins to profile the field. If a units dict is provided, bins are understood to be in the units specified in the dictionary.
Examples
Create a 1d profile. Access bin field from profile.x and field data from profile[<field_name>].
>>> ds = load("DD0046/DD0046") >>> ad = ds.all_data() >>> profile = ad.profile( ... ad, ... [("gas", "density")], ... [("gas", "temperature"), ("gas", "velocity_x")], ... ) >>> print(profile.x) >>> print(profile["gas", "temperature"]) >>> plot = profile.plot()
- ptp(field)¶
Compute the range of values (maximum - minimum) of a field.
This will, in a parallel-aware fashion, compute the “peak-to-peak” of the given field.
- Parameters:
field (string or tuple field name) – The field to average.
- Return type:
Scalar
Examples
>>> rho_range = reg.ptp(("gas", "density"))
- save_as_dataset(filename=None, fields=None)¶
Export a data object to a reloadable yt dataset.
This function will take a data object and output a dataset containing either the fields presently existing or fields given in the
fields
list. The resulting dataset can be reloaded as a yt dataset.- Parameters:
filename (str, optional) – The name of the file to be written. If None, the name will be a combination of the original dataset and the type of data container.
fields (list of string or tuple field names, optional) – If this is supplied, it is the list of fields to be saved to disk. If not supplied, all the fields that have been queried will be saved.
- Returns:
filename – The name of the file that has been created.
- Return type:
Examples
>>> import yt >>> ds = yt.load("enzo_tiny_cosmology/DD0046/DD0046") >>> sp = ds.sphere(ds.domain_center, (10, "Mpc")) >>> fn = sp.save_as_dataset(fields=[("gas", "density"), ("gas", "temperature")]) >>> sphere_ds = yt.load(fn) >>> # the original data container is available as the data attribute >>> print(sds.data["gas", "density"]) [ 4.46237613e-32 4.86830178e-32 4.46335118e-32 ..., 6.43956165e-30 3.57339907e-30 2.83150720e-30] g/cm**3 >>> ad = sphere_ds.all_data() >>> print(ad["gas", "temperature"]) [ 1.00000000e+00 1.00000000e+00 1.00000000e+00 ..., 4.40108359e+04 4.54380547e+04 4.72560117e+04] K
- property selector¶
- set_field_parameter(name, val)¶
Here we set up dictionaries that get passed up and down and ultimately to derived fields.
- std(field, axis=None, weight=None)¶
Compute the standard deviation of a field, optionally along an axis, with a weight.
This will, in a parallel-ware fashion, compute the standard deviation of the given field. If an axis is supplied, it will return a projection, where the weight is also supplied.
By default the weight field will be “ones” or “particle_ones”, depending on the field, resulting in an unweighted standard deviation.
- Parameters:
field (string or tuple field name) – The field to calculate the standard deviation of
axis (string, optional) – If supplied, the axis to compute the standard deviation along (i.e., to project along)
weight (string, optional) – The field to use as a weight.
- Return type:
Scalar or YTProjection.
- sum(field, axis=None)¶
Compute the sum of a field, optionally along an axis.
This will, in a parallel-aware fashion, compute the sum of the given field. If an axis is specified, it will return a projection (using method type “sum”, which does not take into account path length) along that axis.
- Parameters:
field (string or tuple field name) – The field to sum.
axis (string, optional) – If supplied, the axis to sum along.
- Return type:
Either a scalar or a YTProjection.
Examples
>>> total_vol = reg.sum("cell_volume") >>> cell_count = reg.sum(("index", "ones"), axis=("index", "x"))
- property tiles¶
- to_astropy_table(fields)¶
Export region data to a :class:~astropy.table.table.QTable, which is a Table object which is unit-aware. The QTable can then be exported to an ASCII file, FITS file, etc.
See the AstroPy Table docs for more details: http://docs.astropy.org/en/stable/table/
- Parameters:
fields (list of strings or tuple field names) – This is the list of fields to be exported into the QTable.
Examples
>>> sp = ds.sphere("c", (1.0, "Mpc")) >>> t = sp.to_astropy_table([("gas", "density"), ("gas", "temperature")])
- to_dataframe(fields)¶
Export a data object to a
DataFrame
.This function will take a data object and an optional list of fields and export them to a
DataFrame
object. If pandas is not importable, this will raise ImportError.- Parameters:
fields (list of strings or tuple field names) – This is the list of fields to be exported into the DataFrame.
- Returns:
df – The data contained in the object.
- Return type:
Examples
>>> dd = ds.all_data() >>> df = dd.to_dataframe([("gas", "density"), ("gas", "temperature")])
- to_glue(fields, label='yt', data_collection=None)¶
Takes specific fields in the container and exports them to Glue (http://glueviz.org) for interactive analysis. Optionally add a label. If you are already within the Glue environment, you can pass a data_collection object, otherwise Glue will be started.
- volume()¶
Return the volume of the data container. This is found by adding up the volume of the cells with centers in the container, rather than using the geometric shape of the container, so this may vary very slightly from what might be expected from the geometric volume.
- write_out(filename, fields=None, format='%0.16e')¶
Write out the YTDataContainer object in a text file.
This function will take a data object and produce a tab delimited text file containing the fields presently existing and the fields given in the
fields
list.- Parameters:
filename (String) – The name of the file to write to.
fields (List of string, Default = None) – If this is supplied, these fields will be added to the list of fields to be saved to disk. If not supplied, whatever fields presently exist will be used.
format (String, Default = "%0.16e") – Format of numbers to be written in the file.
- Raises:
ValueError – Raised when there is no existing field.
YTException – Raised when field_type of supplied fields is inconsistent with the field_type of existing fields.
Examples
>>> ds = fake_particle_ds() >>> sp = ds.sphere(ds.domain_center, 0.25) >>> sp.write_out("sphere_1.txt") >>> sp.write_out("sphere_2.txt", fields=["cell_volume"])