14. Data Processing and Visualization

This section is intended to explain procedures for carrying out common post-processing tasks with OpenMC. While several utilities of varying complexity are provided to help automate the process, the most powerful capabilities for post-processing derive from use of the Python API.

14.1. Working with State Points

Tally results are saved in both a text file (tallies.out) as well as an HDF5 statepoint file. While the tallies.out file may be fine for simple tallies, in many cases the user requires more information about the tally or the run, or has to deal with a large number of result values (e.g. for mesh tallies). In these cases, extracting data from the statepoint file via the Python API is the preferred method of data analysis and visualization.

14.1.1. Data Extraction

A great deal of information is available in statepoint files (See State Point File Format), all of which is accessible through the Python API. The openmc.StatePoint class can load statepoints and access data as requested; it is used in many of the provided plotting utilities, OpenMC’s regression test suite, and can be used in user-created scripts to carry out manipulations of the data.

An example notebook demonstrates how to extract data from a statepoint using the Python API.

14.1.2. Plotting in 2D

The example notebook also demonstrates how to plot a structured mesh tally in two dimensions using the Python API. One can also use the openmc-plotter application that provides an interactive GUI to explore and plot a much wider variety of tallies.

14.2. Particle Track Visualization

../_images/Tracks.png

OpenMC can dump particle tracks—the position of particles as they are transported through the geometry. There are two ways to make OpenMC output tracks: all particle tracks through a command line argument or specific particle tracks through settings.xml.

Running openmc with the argument -t or --track will cause a track file to be created for every particle transported in the code. Be careful as this will produce as many files as there are source particles in your simulation. To identify a specific particle for which a track should be created, set the Settings.track attribute to a tuple containing the batch, generation, and particle number of the desired particle. For example, to create a track file for particle 4 of batch 1 and generation 2:

settings = openmc.Settings()
settings.track = [(1, 2, 4)]

To specify multiple particles, specify a list of tuples, e.g., if we wanted particles 3 and 4 from batch 1 and generation 2:

settings.track = [(1, 2, 3), (1, 2, 4)]

After running OpenMC (now, without the -t argument), the working directory will contain a file named tracks.h5, which contains a collection of particle tracks. These track files can be converted into VTK poly data files or matplotlib plots with the openmc.Tracks class.

14.3. Source Site Processing

For eigenvalue problems, OpenMC will store information on the fission source sites in the statepoint file by default. For each source site, the weight, position, sampled direction, and sampled energy are stored. To extract this data from a statepoint file, the openmc.statepoint module can be used. An example notebook demontrates how to analyze and plot source information.

14.4. VTK Mesh File Generation

VTK files of OpenMC meshes can be created using the openmc.Mesh.write_data_to_vtk() method. Data can be applied to the elements of the resulting mesh from mesh filter objects. This data can be provided either as a flat array or, in the case of structured meshes (RegularMesh, RectilinearMesh, CylindricalMesh, or SphericalMesh), the data can be shaped with dimensions that match the dimensions of the mesh itself.

OpenMC spherical mesh exported to VTK

For all mesh types, if a flat data array is provided to the mesh, it is expected that the data is ordered in the same ordering as the openmc.Mesh.indices for that mesh object. When providing data directly from a tally, as shown below, a flat array for a given dataset can be passed directly to this method.

# create model above

# create a mesh tally
mesh = openmc.RegularMesh()
mesh.dimension = [10, 20, 30]
mesh.lower_left = [-5, -10, -15]
mesh.upper_right = [5, 10, 15]
mesh_filter = openmc.MeshFilter(mesh)
tally = openmc.Tally()
tally.filters = [mesh_filter]
tally.scores = ['flux']

model.tallies = [tally]
model.run(apply_tally_results=True)

# provide the data as-is to the method
mesh.write_data_to_vtk('flux.vtk', {'flux-mean': tally.mean})

The Tally object also provides a way to expand the dimensions of the mesh filter into a meaningful form where indexing the mesh filter dimensions results in intuitive slicing of structured meshes by setting expand_dims=True when using openmc.Tally.get_reshaped_data(). This reshaping does cause flat indexing of the data to change, however. As noted above, provided datasets are allowed to be shaped so long as such datasets have shapes that match the mesh dimensions. The ability to pass datasets in this way is useful when additional filters are applied to a tally. The example below demonstrates such a case for tally with both a MeshFilter and EnergyFilter applied.

# create model above

# create a mesh tally with energy filter
mesh = openmc.RegularMesh()
mesh.dimension = [10, 20, 30]
mesh.lower_left = [-5, -10, -15]
mesh.upper_right = [5, 10, 15]
mesh_filter = openmc.MeshFilter(mesh)
energy_filter = openmc.EnergyFilter([0.0, 1.0, 20.0e6])
tally = openmc.Tally()
tally.filters = [mesh_filter, energy_filter]
tally.scores = ['flux']

model.tallies = [tally]
model.run(apply_tally_results=True)

# get the data with mesh dimensions expanded, squeeze out length-one dimensions (nuclides, scores)
flux = tally.get_reshaped_data(expand_dims=True).squeeze() # shape: (10, 20, 30, 2)

# write the lowest energy group to a VTK file
mesh.write_data_to_vtk('flux-group1.vtk', datasets={'flux-mean': flux[..., 0]})