yt 2.6 documentation

Loading Generic Array DataΒΆ

[]

Even if your data is not strictly related to fields commonly used in astrophysical codes or your code is not supported yet, you can still feed it to yt to use its advanced visualization and analysis facilities. The only requirement is that your data can be represented as three-dimensional NumPy arrays with a consistent grid structure. What follows are some common examples of loading in generic array data that you may find useful.

Generic Unigrid Data

The simplest case is that of a single grid of data spanning the domain, with one or more fields. The data could be generated from a variety of sources; we'll just give three common examples:

Data generated "on-the-fly"

The most common example is that of data that is generated in memory from the currently running script or notebook.

In [1]:
from yt.imods import *
from yt.utilities.physical_constants import cm_per_kpc, cm_per_mpc

In this example, we'll just create a 3-D array of random floating-point data using NumPy:

In [2]:
arr = np.random.random(size=(64,64,64))

To load this data into yt, we need to assign it a field name, in this case "Density", and place it into a dictionary. Then, we call load_uniform_grid:

In [3]:
data = dict(Density = arr)
bbox = np.array([[-1.5, 1.5], [-1.5, 1.5], [-1.5, 1.5]])
pf = load_uniform_grid(data, arr.shape, cm_per_mpc, bbox=bbox, nprocs=64)

load_uniform_grid takes the following arguments and optional keywords:

  • data : This is a dict of numpy arrays, where the keys are the field names.
  • domain_dimensions : The domain dimensions of the unigrid
  • sim_unit_to_cm : Conversion factor from simulation units to centimeters
  • bbox : Size of computational domain in units sim_unit_to_cm
  • nprocs : If greater than 1, will create this number of subarrays out of data
  • sim_time : The simulation time in seconds
  • periodicity : A tuple of booleans that determines whether the data will be treated as periodic along each axis

This example creates a yt-native parameter file pf that will treat your array as a density field in cubic domain of 3 Mpc edge size (3 * 3.0856e24 cm) and simultaneously divide the domain into nprocs = 64 chunks, so that you can take advantage of the underlying parallelism.

The resulting pf functions exactly like a parameter file from any other dataset--it can be sliced, and we can show the grid boundaries:

In [4]:
slc = SlicePlot(pf, 2, ["Density"])
slc.set_cmap("Density", "Blues")
slc.annotate_grids(cmap=None)
slc.show()

Particle fields are detected as one-dimensional fields. The number of particles is set by the number_of_particles key in data. Particle fields are then added as one-dimensional arrays in a similar manner as the three-dimensional grid fields:

In [5]:
posx_arr = np.random.uniform(low=-1.5, high=1.5, size=10000)
posy_arr = np.random.uniform(low=-1.5, high=1.5, size=10000)
posz_arr = np.random.uniform(low=-1.5, high=1.5, size=10000)
data = dict(Density = np.random.random(size=(64,64,64)), 
            number_of_particles = 10000,
            particle_position_x = posx_arr, 
	        particle_position_y = posy_arr,
	        particle_position_z = posz_arr)
bbox = np.array([[-1.5, 1.5], [-1.5, 1.5], [-1.5, 1.5]])
pf = load_uniform_grid(data, data["Density"].shape, cm_per_mpc, bbox=bbox, nprocs=4)

In this example only the particle position fields have been assigned. number_of_particles must be the same size as the particle arrays. If no particle arrays are supplied then number_of_particles is assumed to be zero. Take a slice, and overlay particle positions:

In [6]:
slc = SlicePlot(pf, "z", ["Density"])
slc.set_cmap("Density", "Blues")
slc.annotate_particles(0.25, p_size=12.0, col="Red")
slc.show()

HDF5 data

HDF5 is a convenient format to store data. If you have unigrid data stored in an HDF5 file, it is possible to load it into memory and then use load_uniform_grid to get it into yt:

In [7]:
import h5py
from yt.config import ytcfg
data_dir = ytcfg.get('yt','test_data_dir')
f = h5py.File(data_dir+"/UnigridData/turb_vels.h5", "r") # Read-only access to the file

The HDF5 file handle's keys correspond to the datasets stored in the file:

In [8]:
print f.keys()
[u'Bx', u'By', u'Bz', u'Density', u'MagneticEnergy', u'Temperature', u'turb_x-velocity', u'turb_y-velocity', u'turb_z-velocity', u'x-velocity', u'y-velocity', u'z-velocity']

We can iterate over the items in the file handle to get the data into a dictionary, which we will then load:

In [9]:
data = {k:v for k,v in f.items()}
bbox = np.array([[-0.5, 0.5], [-0.5, 0.5], [-0.5, 0.5]])
pf = load_uniform_grid(data, data["Density"].shape, 250.*cm_per_kpc, bbox=bbox, nprocs=8, periodicity=(False,False,False))

In this case, the data came from a simulation which was 250 kpc on a side. An example projection of two fields:

In [10]:
prj = ProjectionPlot(pf, "z", ["z-velocity","Temperature"], weight_field="Density")
prj.show()