1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248
|
.. _saving_data:
Saving Reloadable Data
======================
Most of the data loaded into or generated with yt can be saved to a
format that can be reloaded as a first-class dataset. This includes
the following:
* geometric data containers (regions, spheres, disks, rays, etc.)
* grid data containers (covering grids, arbitrary grids, fixed
resolution buffers)
* spatial plots (projections, slices, cutting planes)
* profiles
* generic array data
In the case of projections, slices, and profiles, reloaded data can be
used to remake plots. For information on this, see :ref:`remaking-plots`.
.. _saving-data-containers:
Geometric Data Containers
-------------------------
Data from geometric data containers can be saved with the
:func:`~yt.data_objects.data_containers.YTDataContainer.save_as_dataset` function.
.. notebook-cell::
import yt
ds = yt.load("enzo_tiny_cosmology/DD0046/DD0046")
sphere = ds.sphere([0.5] * 3, (10, "Mpc"))
fn = sphere.save_as_dataset(fields=[("gas", "density"), ("all", "particle_mass")])
print(fn)
This function will return the name of the file to which the dataset
was saved. The filename will be a combination of the name of the
original dataset and the type of data container. Optionally, a
specific filename can be given with the ``filename`` keyword. If no
fields are given, the fields that have previously been queried will
be saved.
The newly created dataset can be loaded like all other supported
data through ``yt.load``. Once loaded, field data can be accessed
through the traditional data containers or through the ``data``
attribute, which will be a data container configured like the
original data container used to make the dataset. Grid data is
accessed by the ``grid`` data type and particle data is accessed
with the original particle type. As with the original dataset, grid
positions and cell sizes are accessible with, for example,
("grid", "x") and ("grid", "dx"). Particle positions are
accessible as (<particle_type>, "particle_position_x"). All original
simulation parameters are accessible in the ``parameters``
dictionary, normally associated with all datasets.
.. code-block:: python
sphere_ds = yt.load("DD0046_sphere.h5")
# use the original data container
print(sphere_ds.data["grid", "density"])
# create a new data container
ad = sphere_ds.all_data()
# grid data
print(ad["grid", "density"])
print(ad["grid", "x"])
print(ad["grid", "dx"])
# particle data
print(ad["all", "particle_mass"])
print(ad["all", "particle_position_x"])
Note that because field data queried from geometric containers is
returned as unordered 1D arrays, data container datasets are treated,
effectively, as particle data. Thus, 3D indexing of grid data from
these datasets is not possible.
.. _saving-grid-data-containers:
Grid Data Containers
--------------------
Data containers that return field data as multidimensional arrays
can be saved so as to preserve this type of access. This includes
covering grids, arbitrary grids, and fixed resolution buffers.
Saving data from these containers works just as with geometric data
containers. Field data can be accessed through geometric data
containers.
.. code-block:: python
cg = ds.covering_grid(level=0, left_edge=[0.25] * 3, dims=[16] * 3)
fn = cg.save_as_dataset(fields=[("gas", "density"), ("all", "particle_mass")])
cg_ds = yt.load(fn)
ad = cg_ds.all_data()
print(ad["grid", "density"])
Multidimensional indexing of field data is also available through
the ``data`` attribute.
.. code-block:: python
print(cg_ds.data["grid", "density"])
Fixed resolution buffers work just the same.
.. code-block:: python
my_proj = ds.proj(("gas", "density"), "x", weight_field=("gas", "density"))
frb = my_proj.to_frb(1.0, (800, 800))
fn = frb.save_as_dataset(fields=[("gas", "density")])
frb_ds = yt.load(fn)
print(frb_ds.data["gas", "density"])
.. _saving-spatial-plots:
Spatial Plots
-------------
Spatial plots, such as projections, slices, and off-axis slices
(cutting planes) can also be saved and reloaded.
.. code-block:: python
proj = ds.proj(("gas", "density"), "x", weight_field=("gas", "density"))
proj.save_as_dataset()
Once reloaded, they can be handed to their associated plotting
functions to make images.
.. code-block:: python
proj_ds = yt.load("DD0046_proj.h5")
p = yt.ProjectionPlot(proj_ds, "x", ("gas", "density"), weight_field=("gas", "density"))
p.save()
.. _saving-profile-data:
Profiles
--------
Profiles created with :func:`~yt.data_objects.profiles.create_profile`,
:class:`~yt.visualization.profile_plotter.ProfilePlot`, and
:class:`~yt.visualization.profile_plotter.PhasePlot` can be saved with
the :func:`~yt.data_objects.profiles.save_as_dataset` function, which
works just as above. Profile datasets are a type of non-spatial grid
datasets. Geometric selection is not possible, but data can be
accessed through the ``.data`` attribute.
.. notebook-cell::
import yt
ds = yt.load("enzo_tiny_cosmology/DD0046/DD0046")
ad = ds.all_data()
profile_2d = yt.create_profile(ad, [("gas", "density"), ("gas", "temperature")],
("gas", "mass"), weight_field=None,
n_bins=(128, 128))
profile_2d.save_as_dataset()
prof_2d_ds = yt.load("DD0046_Profile2D.h5")
print (prof_2d_ds.data["gas", "mass"])
The x, y (if at least 2D), and z (if 3D) bin fields can be accessed as 1D
arrays with "x", "y", and "z".
.. code-block:: python
print(prof_2d_ds.data["gas", "x"])
The bin fields can also be returned with the same shape as the profile
data by accessing them with their original names. This allows for
boolean masking of profile data using the bin fields.
.. code-block:: python
# density is the x bin field
print(prof_2d_ds.data["gas", "density"])
For 1, 2, and 3D profile datasets, a fake profile object will be
constructed by accessing the ".profile" attribute. This is used
primarily in the case of 1 and 2D profiles to create figures using
:class:`~yt.visualization.profile_plotter.ProfilePlot` and
:class:`~yt.visualization.profile_plotter.PhasePlot`.
.. code-block:: python
p = yt.PhasePlot(
prof_2d_ds.data,
("gas", "density"),
("gas", "temperature"),
("gas", "mass"),
weight_field=None,
)
p.save()
.. _saving-array-data:
Generic Array Data
------------------
Generic arrays can be saved and reloaded as non-spatial data using
the :func:`~yt.frontends.ytdata.utilities.save_as_dataset` function,
also available as ``yt.save_as_dataset``. As with profiles, geometric
selection is not possible, but the data can be accessed through the
``.data`` attribute.
.. notebook-cell::
import yt
ds = yt.load("enzo_tiny_cosmology/DD0046/DD0046")
region = ds.box([0.25]*3, [0.75]*3)
sphere = ds.sphere(ds.domain_center, (10, "Mpc"))
my_data = {}
my_data["region_density"] = region["gas", "density"]
my_data["sphere_density"] = sphere["gas", "density"]
yt.save_as_dataset(ds, "test_data.h5", my_data)
array_ds = yt.load("test_data.h5")
print (array_ds.data["data", "region_density"])
print (array_ds.data["data", "sphere_density"])
Array data can be saved with or without a dataset loaded. If no
dataset has been loaded, as fake dataset can be provided as a
dictionary.
.. notebook-cell::
import numpy as np
import yt
my_data = {"density": yt.YTArray(np.random.random(10), "g/cm**3"),
"temperature": yt.YTArray(np.random.random(10), "K")}
fake_ds = {"current_time": yt.YTQuantity(10, "Myr")}
yt.save_as_dataset(fake_ds, "random_data.h5", my_data)
new_ds = yt.load("random_data.h5")
print (new_ds.data["data", "density"])
|