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.
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:
The most common example is that of data that is generated in memory from the currently running script or notebook.
In [ ]:
import yt
import numpy as np
In this example, we'll just create a 3-D array of random floating-point data using NumPy:
In [ ]:
arr = np.random.random(size=(64,64,64))
To load this data into yt, we need associate it with a field. The data dictionary consists of one or more fields, each consisting of a tuple of a NumPy array and a unit string. Then, we can call load_uniform_grid:
In [ ]:
data = dict(density = (arr, "g/cm**3"))
bbox = np.array([[-1.5, 1.5], [-1.5, 1.5], [-1.5, 1.5]])
ds = yt.load_uniform_grid(data, arr.shape, length_unit="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 namesdomain_dimensions : The domain dimensions of the unigridlength_unit : The unit that corresponds to code_length, can be a string, tuple, or floating-point numberbbox : Size of computational domain in units of code_lengthnprocs : If greater than 1, will create this number of subarrays out of datasim_time : The simulation time in secondsmass_unit : The unit that corresponds to code_mass, can be a string, tuple, or floating-point numbertime_unit : The unit that corresponds to code_time, can be a string, tuple, or floating-point numbervelocity_unit : The unit that corresponds to code_velocitymagnetic_unit : The unit that corresponds to code_magnetic, i.e. the internal units used to represent magnetic field strengths.periodicity : A tuple of booleans that determines whether the data will be treated as periodic along each axisThis example creates a yt-native dataset ds that will treat your array as a
density field in cubic domain of 3 Mpc edge size and simultaneously divide the
domain into nprocs = 64 chunks, so that you can take advantage
of the underlying parallelism.
The optional unit keyword arguments allow for the default units of the dataset to be set. They can be:
length_unit="Mpc"mass_unit=(1.0e14, "Msun")time_unit=3.1557e13In the latter case, the unit is assumed to be cgs.
The resulting ds functions exactly like a dataset like any other yt can handle--it can be sliced, and we can show the grid boundaries:
In [ ]:
slc = yt.SlicePlot(ds, "z", ["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 [ ]:
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)), "Msun/kpc**3"),
number_of_particles = 10000,
particle_position_x = (posx_arr, 'code_length'),
particle_position_y = (posy_arr, 'code_length'),
particle_position_z = (posz_arr, 'code_length'))
bbox = np.array([[-1.5, 1.5], [-1.5, 1.5], [-1.5, 1.5]])
ds = yt.load_uniform_grid(data, data["density"][0].shape, length_unit=(1.0, "Mpc"), mass_unit=(1.0,"Msun"),
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 [ ]:
slc = yt.SlicePlot(ds, "z", ["density"])
slc.set_cmap("density", "Blues")
slc.annotate_particles(0.25, p_size=12.0, col="Red")
slc.show()
In a similar fashion to unigrid data, data gridded into rectangular patches at varying levels of resolution may also be loaded into yt. In this case, a list of grid dictionaries should be provided, with the requisite information about each grid's properties. This example sets up two grids: a top-level grid (level == 0) covering the entire domain and a subgrid at level == 1.
In [ ]:
grid_data = [
dict(left_edge = [0.0, 0.0, 0.0],
right_edge = [1.0, 1.0, 1.0],
level = 0,
dimensions = [32, 32, 32]),
dict(left_edge = [0.25, 0.25, 0.25],
right_edge = [0.75, 0.75, 0.75],
level = 1,
dimensions = [32, 32, 32])
]
We'll just fill each grid with random density data, with a scaling with the grid refinement level.
In [ ]:
for g in grid_data:
g["density"] = (np.random.random(g["dimensions"]) * 2**g["level"], "g/cm**3")
Particle fields are supported by adding 1-dimensional arrays to each grid and
setting the number_of_particles key in each grid's dict. If a grid has no particles, set number_of_particles = 0, but the particle fields still have to be defined since they are defined elsewhere; set them to empty NumPy arrays:
In [ ]:
grid_data[0]["number_of_particles"] = 0 # Set no particles in the top-level grid
grid_data[0]["particle_position_x"] = (np.array([]), "code_length") # No particles, so set empty arrays
grid_data[0]["particle_position_y"] = (np.array([]), "code_length")
grid_data[0]["particle_position_z"] = (np.array([]), "code_length")
grid_data[1]["number_of_particles"] = 1000
grid_data[1]["particle_position_x"] = (np.random.uniform(low=0.25, high=0.75, size=1000), "code_length")
grid_data[1]["particle_position_y"] = (np.random.uniform(low=0.25, high=0.75, size=1000), "code_length")
grid_data[1]["particle_position_z"] = (np.random.uniform(low=0.25, high=0.75, size=1000), "code_length")
Then, call load_amr_grids:
In [ ]:
ds = yt.load_amr_grids(grid_data, [32, 32, 32])
load_amr_grids also takes the same keywords bbox and sim_time as load_uniform_grid. We could have also specified the length, time, velocity, and mass units in the same manner as before. Let's take a slice:
In [ ]:
slc = yt.SlicePlot(ds, "z", ["density"])
slc.annotate_particles(0.25, p_size=15.0, col="Pink")
slc.show()
load_amr_grids assumes that particle positions associated with one grid are not bounded within another grid at a higher level, so this must be ensured by the user prior to loading the grid data.