1.. _generating-processed-data: 2 3Generating Processed Data 4========================= 5 6Although yt provides a number of built-in visualization methods that can 7process data and construct from that plots, it is often useful to generate the 8data by hand and construct plots which can then be combined with other plots, 9modified in some way, or even (gasp) created and modified in some other tool or 10program. 11 12.. _exporting-container-data: 13 14Exporting Container Data 15------------------------ 16 17Fields from data containers such as regions, spheres, cylinders, etc. can be exported 18tabular format using either a :class:`~pandas.DataFrame` or an :class:`~astropy.table.QTable`. 19 20To export to a :class:`~pandas.DataFrame`, use 21:meth:`~yt.data_objects.data_containers.YTDataContainer.to_dataframe`: 22 23.. code-block:: python 24 25 sp = ds.sphere("c", (0.2, "unitary")) 26 df2 = sp.to_dataframe([("gas", "density"), ("gas", "temperature")]) 27 28To export to a :class:`~astropy.table.QTable`, use 29:meth:`~yt.data_objects.data_containers.YTDataContainer.to_astropy_table`: 30 31.. code-block:: python 32 33 sp = ds.sphere("c", (0.2, "unitary")) 34 at2 = sp.to_astropy_table(fields=[("gas", "density"), ("gas", "temperature")]) 35 36For exports to :class:`~pandas.DataFrame` objects, the unit information is lost, but for 37exports to :class:`~astropy.table.QTable` objects, the :class:`~yt.units.yt_array.YTArray` 38objects are converted to :class:`~astropy.units.Quantity` objects. 39 40.. _generating-2d-image-arrays: 41 422D Image Arrays 43--------------- 44 45When making a slice, a projection or an oblique slice in yt, the resultant 46:class:`~yt.data_objects.data_containers.YTSelectionContainer2D` object is created and 47contains flattened arrays of the finest available data. This means a set of 48arrays for the x, y, (possibly z), dx, dy, (possibly dz) and data values, for 49every point that constitutes the object. 50 51 52This presents something of a challenge for visualization, as it will require 53the transformation of a variable mesh of points consisting of positions and 54sizes into a fixed-size array that appears like an image. This process is that 55of pixelization, which yt handles transparently internally. You can access 56this functionality by constructing a 57:class:`~yt.visualization.fixed_resolution.FixedResolutionBuffer` and supplying 58to it your :class:`~yt.data_objects.data_containers.YTSelectionContainer2D` 59object, as well as some information about how you want the final image to look. 60You can specify both the bounds of the image (in the appropriate x-y plane) and 61the resolution of the output image. You can then have yt pixelize any field 62you like. 63 64.. note:: In previous versions of yt, there was a special class of 65 FixedResolutionBuffer for off-axis slices. This is no longer 66 necessary. 67 68To create :class:`~yt.data_objects.data_containers.YTSelectionContainer2D` objects, you can 69access them as described in :ref:`data-objects`, specifically the section 70:ref:`available-objects`. Here is an example of how to window into a slice 71of resolution(512, 512) with bounds of (0.3, 0.5) and (0.6, 0.8). The next 72step is to generate the actual 2D image array, which is accomplished by 73accessing the desired field. 74 75.. code-block:: python 76 77 sl = ds.slice(0, 0.5) 78 frb = FixedResolutionBuffer(sl, (0.3, 0.5, 0.6, 0.8), (512, 512)) 79 my_image = frb["density"] 80 81This image may then be used in a hand-constructed Matplotlib image, for instance using 82:func:`~matplotlib.pyplot.imshow`. 83 84The buffer arrays can be saved out to disk in either HDF5 or FITS format: 85 86.. code-block:: python 87 88 frb.save_as_dataset("my_images.h5", fields=[("gas", "density"), ("gas", "temperature")]) 89 frb.export_fits( 90 "my_images.fits", 91 fields=[("gas", "density"), ("gas", "temperature")], 92 clobber=True, 93 units="kpc", 94 ) 95 96In the HDF5 case, the created file can be reloaded just like a regular dataset with 97``yt.load`` and will, itself, be a first-class dataset. For more information on 98this, see :ref:`saving-grid-data-containers`. 99In the FITS case, there is an option for setting the ``units`` of the coordinate system in 100the file. If you want to overwrite a file with the same name, set ``clobber=True``. 101 102The :class:`~yt.visualization.fixed_resolution.FixedResolutionBuffer` can even be exported 103as a 2D dataset itself, which may be operated on in the same way as any other dataset in yt: 104 105.. code-block:: python 106 107 ds_frb = frb.export_dataset( 108 fields=[("gas", "density"), ("gas", "temperature")], nprocs=8 109 ) 110 sp = ds_frb.sphere("c", (100.0, "kpc")) 111 112where the ``nprocs`` parameter can be used to decompose the image into ``nprocs`` number of grids. 113 114.. _generating-profiles-and-histograms: 115 116Profiles and Histograms 117----------------------- 118 119Profiles and histograms can also be generated using the 120:class:`~yt.visualization.profile_plotter.ProfilePlot` and 121:class:`~yt.visualization.profile_plotter.PhasePlot` functions 122(described in :ref:`how-to-make-1d-profiles` and 123:ref:`how-to-make-2d-profiles`). These generate profiles transparently, but the 124objects they handle and create can be handled manually, as well, for more 125control and access. The :func:`~yt.data_objects.profiles.create_profile` function 126can be used to generate 1, 2, and 3D profiles. 127 128Profile objects can be created from any data object (see :ref:`data-objects`, 129specifically the section :ref:`available-objects` for more information) and are 130best thought of as distribution calculations. They can either sum up or average 131one quantity with respect to one or more other quantities, and they do this over 132all the data contained in their source object. When calculating average values, 133the standard deviation will also be calculated. 134 135To generate a profile, one need only specify the binning fields and the field 136to be profiled. The binning fields are given together in a list. The 137:func:`~yt.data_objects.profiles.create_profile` function will guess the 138dimensionality of the profile based on the number of fields given. For example, 139a one-dimensional profile of the mass-weighted average temperature as a function of 140density within a sphere can be created in the following way: 141 142.. code-block:: python 143 144 import yt 145 146 ds = yt.load("galaxy0030/galaxy0030") 147 source = ds.sphere("c", (10, "kpc")) 148 profile = source.profile( 149 [("gas", "density")], # the bin field 150 [ 151 ("gas", "temperature"), # profile field 152 ("gas", "radial_velocity"), 153 ], # profile field 154 weight_field=("gas", "mass"), 155 ) 156 157The binning, weight, and profile data can now be access as: 158 159.. code-block:: python 160 161 print(profile.x) # bin field 162 print(profile.weight) # weight field 163 print(profile["gas", "temperature"]) # profile field 164 print(profile["gas", "radial_velocity"]) # profile field 165 166The ``profile.used`` attribute gives a boolean array of the bins which actually 167have data. 168 169.. code-block:: python 170 171 print(profile.used) 172 173If a weight field was given, the profile data will represent the weighted mean 174of a field. In this case, the weighted standard deviation will be calculated 175automatically and can be access via the ``profile.standard_deviation`` 176attribute. 177 178.. code-block:: python 179 180 print(profile.standard_deviation["gas", "temperature"]) 181 182A two-dimensional profile of the total gas mass in bins of density and 183temperature can be created as follows: 184 185.. code-block:: python 186 187 profile2d = source.profile( 188 [ 189 ("gas", "density"), 190 ("gas", "temperature"), 191 ], # the x bin field # the y bin field 192 [("gas", "mass")], # the profile field 193 weight_field=None, 194 ) 195 196Accessing the x, y, and profile fields work just as with one-dimensional profiles: 197 198.. code-block:: python 199 200 print(profile2d.x) 201 print(profile2d.y) 202 print(profile2d["gas", "mass"]) 203 204One of the more interesting things that is enabled with this approach is 205the generation of 1D profiles that correspond to 2D profiles. For instance, a 206phase plot that shows the distribution of mass in the density-temperature 207plane, with the average temperature overplotted. The 208:func:`~matplotlib.pyplot.pcolormesh` function can be used to manually plot 209the 2D profile. If you want to generate a default profile plot, you can simply 210call::: 211 212 profile.plot() 213 214Three-dimensional profiles can be generated and accessed following 215the same procedures. Additional keyword arguments are available to control 216the following for each of the bin fields: the number of bins, min and max, units, 217whether to use a log or linear scale, and whether or not to do accumulation to 218create a cumulative distribution function. For more information, see the API 219documentation on the :func:`~yt.data_objects.profiles.create_profile` function. 220 221For custom bins the other keyword arguments can be overriden using the 222``override_bins`` keyword argument. This accepts a dictionary with an array 223for each bin field or ``None`` to use the default settings. 224 225.. code-block:: python 226 227 custom_bins = np.array([1e-27, 1e-25, 2e-25, 5e-25, 1e-23]) 228 profile2d = source.profile( 229 [("gas", "density"), ("gas", "temperature")], 230 [("gas", "mass")], 231 override_bins={("gas", "density"): custom_bins, ("gas", "temperature"): None}, 232 ) 233 234.. _profile-dataframe-export: 235 236Exporting Profiles to DataFrame 237^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ 238 239One-dimensional profile data can be exported to a :class:`~pandas.DataFrame` object 240using the :meth:`yt.data_objects.profiles.Profile1D.to_dataframe` method. Bins which 241do not have data will have their fields filled with ``NaN``, except for the bin field 242itself. If you only want to export the bins which are used, set ``only_used=True``. 243 244.. code-block:: python 245 246 # Adds all of the data to the DataFrame, but non-used bins are filled with NaNs 247 df = profile.to_dataframe() 248 # Only adds the used bins to the DataFrame 249 df_used = profile.to_dataframe(only_used=True) 250 # Only adds the density and temperature fields 251 df2 = profile.to_dataframe(fields=[("gas", "density"), ("gas", "temperature")]) 252 253The :class:`~pandas.DataFrame` can then analyzed and/or written to disk using pandas 254methods. Note that unit information is lost in this export. 255 256.. _profile-astropy-export: 257 258Exporting Profiles to QTable 259^^^^^^^^^^^^^^^^^^^^^^^^^^^^ 260 261One-dimensional profile data also can be exported to an AstroPy :class:`~astropy.table.QTable` 262object. This table can then be written to disk in a number of formats, such as ASCII text 263or FITS files, and manipulated in a number of ways. Bins which do not have data 264will have their mask values set to ``False``. If you only want to export the bins 265which are used, set ``only_used=True``. Units are preserved in the table by converting 266each :class:`~yt.units.yt_array.YTArray` to an :class:`~astropy.units.Quantity`. 267 268To export the 1D profile to a Table object, simply call 269:meth:`yt.data_objects.profiles.Profile1D.to_astropy_table`: 270 271.. code-block:: python 272 273 # Adds all of the data to the Table, but non-used bins are masked 274 t = profile.to_astropy_table() 275 # Only adds the used bins to the Table 276 t_used = profile.to_astropy_table(only_used=True) 277 # Only adds the density and temperature fields 278 t2 = profile.to_astropy_table(fields=[("gas", "density"), ("gas", "temperature")]) 279 280.. _generating-line-queries: 281 282Line Queries and Planar Integrals 283--------------------------------- 284 285To calculate the values along a line connecting two points in a simulation, you 286can use the object :class:`~yt.data_objects.selection_data_containers.YTRay`, 287accessible as the ``ray`` property on a index. (See :ref:`data-objects` 288for more information on this.) To do so, you can supply two points and access 289fields within the returned object. For instance, this code will generate a ray 290between the points (0.3, 0.5, 0.9) and (0.1, 0.8, 0.5) and examine the density 291along that ray: 292 293.. code-block:: python 294 295 ray = ds.ray((0.3, 0.5, 0.9), (0.1, 0.8, 0.5)) 296 print(ray["gas", "density"]) 297 298The points are not ordered, so you may need to sort the data (see the 299example in the 300:class:`~yt.data_objects.selection_data_containers.YTRay` docs). Also 301note, the ray is traversing cells of varying length, as well as 302taking a varying distance to cross each cell. To determine the 303distance traveled by the ray within each cell (for instance, for 304integration) the field ``dt`` is available; this field will sum to 3051.0, as the ray's path will be normalized to 1.0, independent of how 306far it travels through the domain. To determine the value of ``t`` at 307which the ray enters each cell, the field ``t`` is available. For 308instance: 309 310.. code-block:: python 311 312 print(ray["dts"].sum()) 313 print(ray["t"]) 314 315These can be used as inputs to, for instance, the Matplotlib function 316:func:`~matplotlib.pyplot.plot`, or they can be saved to disk. 317 318The volume rendering functionality in yt can also be used to calculate 319off-axis plane integrals, using the 320:class:`~yt.visualization.volume_rendering.transfer_functions.ProjectionTransferFunction` 321in a manner similar to that described in :ref:`volume_rendering`. 322 323.. _generating-xarray: 324 325Regular Grids to xarray 326----------------------- 327 328Objects that subclass from 329:class:`~yt.data_objects.construction_data_containers.YTCoveringGrid` are able 330to export to `xarray <https://xarray.pydata.org/>`_. This enables 331interoperability with anything that can take xarray data. The classes that can do this are 332:class:`~yt.data_objects.construction_data_containers.YTCoveringGrid`, 333:class:`~yt.data_objects.construction_data_containers.YTArbitraryGrid`, and 334:class:`~yt.data_objects.construction_data_containers.YTSmoothedCoveringGrid`. For example, you can: 335 336.. code-block:: python 337 338 grid = ds.r[::256j, ::256j, ::256j] 339 obj = grid.to_xarray(fields=[("gas", "density"), ("gas", "temperature")]) 340 341The returned object, ``obj``, will now have the correct labelled axes and so forth. 342