openmc.deplete.Results

class openmc.deplete.Results(filename='depletion_results.h5')[source]

Results from a depletion simulation

The Results class acts as a list that stores the results from each depletion step and provides extra methods for interrogating these results.

Changed in version 0.13.1: Name changed from ResultsList to Results

Parameters

filename (str, optional) – Path to depletion result file

export_to_materials(burnup_index: int, nuc_with_data: Optional[Iterable[str]] = None, path: Union[str, PathLike] = 'materials.xml') Materials[source]

Return openmc.Materials object based on results at a given step

New in version 0.12.1.

Parameters
  • burn_index (int) – Index of burnup step to evaluate. See also: get_step_where for obtaining burnup step indices from other data such as the time.

  • nuc_with_data (Iterable of str, optional) – Nuclides to include in resulting materials. This can be specified if not all nuclides appearing in depletion results have associated neutron cross sections, and as such cannot be used in subsequent transport calculations. If not provided, nuclides from the cross_sections element of materials.xml will be used. If that element is not present, nuclides from openmc.config[‘cross_sections’] will be used.

  • path (PathLike) –

    Path to materials XML file to read. Defaults to ‘materials.xml’.

    New in version 0.13.3.

Returns

mat_file – A modified Materials instance containing depleted material data and original isotopic compositions of non-depletable materials

Return type

Materials

classmethod from_hdf5(filename: Union[str, PathLike])[source]

Load in depletion results from a previous file

Parameters

filename (str) – Path to depletion result file

Returns

New instance of depletion results

Return type

Results

get_activity(mat: Union[Material, str], units: str = 'Bq/cm3', by_nuclide: bool = False, volume: Optional[float] = None) Tuple[ndarray, Union[ndarray, List[dict]]][source]

Get activity of material over time.

New in version 0.14.0.

Parameters
  • mat (openmc.Material, str) – Material object or material id to evaluate

  • units ({'Bq', 'Bq/g', 'Bq/cm3'}) – Specifies the type of activity to return, options include total activity [Bq], specific [Bq/g] or volumetric activity [Bq/cm3].

  • by_nuclide (bool) – Specifies if the activity should be returned for the material as a whole or per nuclide. Default is False.

  • volume (float, optional) – Volume of the material. If not passed, defaults to using the Material.volume attribute.

Returns

  • times (numpy.ndarray) – Array of times in [s]

  • activities (numpy.ndarray or List[dict]) – Array of total activities if by_nuclide = False (default) or list of dictionaries of activities by nuclide if by_nuclide = True.

get_atoms(mat: Union[Material, str], nuc: str, nuc_units: str = 'atoms', time_units: str = 's') Tuple[ndarray, ndarray][source]

Get number of nuclides over time from a single material

Parameters
  • mat (openmc.Material, str) – Material object or material id to evaluate

  • nuc (str) – Nuclide name to evaluate

  • nuc_units ({"atoms", "atom/b-cm", "atom/cm3"}, optional) –

    Units for the returned concentration. Default is "atoms"

    New in version 0.12.

  • time_units ({"s", "min", "h", "d", "a"}, optional) –

    Units for the returned time array. Default is "s" to return the value in seconds. Other options are minutes "min", hours "h", days "d", and Julian years "a".

    New in version 0.12.

Returns

  • times (numpy.ndarray) – Array of times in units of time_units

  • concentrations (numpy.ndarray) – Concentration of specified nuclide in units of nuc_units

get_decay_heat(mat: Union[Material, str], units: str = 'W', by_nuclide: bool = False, volume: Optional[float] = None) Tuple[ndarray, Union[ndarray, List[dict]]][source]

Get decay heat of material over time.

New in version 0.14.0.

Parameters
  • mat (openmc.Material, str) – Material object or material id to evaluate.

  • units ({'W', 'W/g', 'W/cm3'}) – Specifies the units of decay heat to return. Options include total heat [W], specific [W/g] or volumetric heat [W/cm3].

  • by_nuclide (bool) – Specifies if the decay heat should be returned for the material as a whole or per nuclide. Default is False.

  • volume (float, optional) – Volume of the material. If not passed, defaults to using the Material.volume attribute.

Returns

  • times (numpy.ndarray) – Array of times in [s]

  • decay_heat (numpy.ndarray or List[dict]) – Array of total decay heat values if by_nuclide = False (default) or list of dictionaries of decay heat values by nuclide if by_nuclide = True.

get_depletion_time() ndarray[source]

Return an array of the average time to deplete a material

Note

The return value will have one fewer values than several other methods, such as get_keff(), because no depletion is performed at the final transport stage.

Returns

times – Vector of average time to deplete a single material across all processes and materials.

Return type

numpy.ndarray

get_keff(time_units: str = 's') Tuple[ndarray, ndarray][source]

Evaluates the eigenvalue from a results list.

New in version 0.13.1.

Parameters

time_units ({"s", "d", "min", "h", "a"}, optional) – Desired units for the times array. Options are seconds "s", minutes "min", hours "h", days "d", and Julian years "a".

Returns

  • times (numpy.ndarray) – Array of times in specified units

  • eigenvalues (numpy.ndarray) – k-eigenvalue at each time. Column 0 contains the eigenvalue, while column 1 contains the associated uncertainty

get_mass(mat: Union[Material, str], nuc: str, mass_units: str = 'g', time_units: str = 's') Tuple[ndarray, ndarray][source]

Get mass of nuclides over time from a single material

New in version 0.14.0.

Parameters
  • mat (openmc.Material, str) – Material object or material id to evaluate

  • nuc (str) – Nuclide name to evaluate

  • mass_units ({"g", "g/cm3", "kg"}, optional) – Units for the returned mass.

  • time_units ({"s", "min", "h", "d", "a"}, optional) – Units for the returned time array. Default is "s" to return the value in seconds. Other options are minutes "min", hours "h", days "d", and Julian years "a".

Returns

  • times (numpy.ndarray) – Array of times in units of time_units

  • mass (numpy.ndarray) – Mass of specified nuclide in units of mass_units

get_reaction_rate(mat: Union[Material, str], nuc: str, rx: str) Tuple[ndarray, ndarray][source]

Get reaction rate in a single material/nuclide over time

Parameters
  • mat (openmc.Material, str) – Material object or material id to evaluate

  • nuc (str) – Nuclide name to evaluate

  • rx (str) – Reaction rate to evaluate

Returns

  • times (numpy.ndarray) – Array of times in [s]

  • rates (numpy.ndarray) – Array of reaction rates

get_step_where(time, time_units: str = 'd', atol: float = 1e-06, rtol: float = 0.001) int[source]

Return the index closest to a given point in time

In the event time lies exactly between two points, the lower index will be returned. It is possible that the index will be at most one past the point in time requested, but only according to tolerances requested.

Passing atol=math.inf and rtol=math.inf will return the closest index to the requested point.

New in version 0.12.1.

Parameters
  • time (float) – Desired point in time

  • time_units ({"s", "d", "min", "h", "a"}, optional) – Units on time. Default: days "d". Other options are seconds "s", minutes "min", hours "h" and Julian years "a".

  • atol (float, optional) – Absolute tolerance (in time_units) if time is not found.

  • rtol (float, optional) – Relative tolerance if time is not found.

Return type

int

get_times(time_units: str = 'd') ndarray[source]

Return the points in time that define the depletion schedule

New in version 0.12.1.

Parameters

time_units ({"s", "d", "min", "h", "a"}, optional) – Return the vector in these units. Default is to convert to days "d". Other options are seconds "s", minutes "min", hours "h", days "d", and Julian years "a".

Returns

1-D vector of time points

Return type

numpy.ndarray