openmc.deplete.Results

class openmc.deplete.Results[source]

Output of a depletion run

Variables:
  • k (list of (float, float)) – Eigenvalue and uncertainty for each substep.
  • time (list of float) – Time at beginning, end of step, in seconds.
  • source_rate (float) – Source rate during timestep in [W] or [neutron/sec]
  • n_mat (int) – Number of mats.
  • n_nuc (int) – Number of nuclides.
  • rates (list of ReactionRates) – The reaction rates for each substep.
  • volume (OrderedDict of str to float) – Dictionary mapping mat id to volume.
  • mat_to_ind (OrderedDict of str to int) – A dictionary mapping mat ID as string to index.
  • nuc_to_ind (OrderedDict of str to int) – A dictionary mapping nuclide name as string to index.
  • mat_to_hdf5_ind (OrderedDict of str to int) – A dictionary mapping mat ID as string to global index.
  • n_hdf5_mats (int) – Number of materials in entire geometry.
  • n_stages (int) – Number of stages in simulation.
  • data (numpy.ndarray) – Atom quantity, stored by stage, mat, then by nuclide.
  • proc_time (int) – Average time spent depleting a material across all materials and processes
allocate(volume, nuc_list, burn_list, full_burn_list, stages)[source]

Allocates memory of Results.

Parameters:
  • volume (dict of str float) – Volumes corresponding to materials in full_burn_dict
  • nuc_list (list of str) – A list of all nuclide names. Used for sorting the simulation.
  • burn_list (list of int) – A list of all mat IDs to be burned. Used for sorting the simulation.
  • full_burn_list (list of str) – List of all burnable material IDs
  • stages (int) – Number of stages in simulation.
distribute(local_materials, ranges)[source]

Create a new object containing data for distributed materials

Parameters:
  • local_materials (iterable of str) – Materials for this process
  • ranges (iterable of int) – Slice-like object indicating indicies of local_materials in the material dimension of data and each element in rates
Returns:

New results object

Return type:

Results

export_to_hdf5(filename, step)[source]

Export results to an HDF5 file

Parameters:
  • filename (str) – The filename to write to
  • step (int) – What step is this?
classmethod from_hdf5(handle, step)[source]

Loads results object from HDF5.

Parameters:
  • handle (h5py.File or h5py.Group) – An HDF5 file or group type to load from.
  • step (int) – Index for depletion step
static save(op, x, op_results, t, source_rate, step_ind, proc_time=None)[source]

Creates and writes depletion results to disk

Parameters:
  • op (openmc.deplete.TransportOperator) – The operator used to generate these results.
  • x (list of list of numpy.array) – The prior x vectors. Indexed [i][cell] using the above equation.
  • op_results (list of openmc.deplete.OperatorResult) – Results of applying transport operator
  • t (list of float) – Time indices.
  • source_rate (float) – Source rate during time step in [W] or [neutron/sec]
  • step_ind (int) – Step index.
  • proc_time (float or None) – Total process time spent depleting materials. This may be process-dependent and will be reduced across MPI processes.
transfer_volumes(geometry)[source]

Transfers volumes from depletion results to geometry

Parameters:geometry (OpenMC geometry to be used in a depletion restart) – calculation