openmc.deplete.ResultsList¶
-
class
openmc.deplete.
ResultsList
[source]¶ A list of openmc.deplete.Results objects
It is recommended to use
from_hdf5()
over direct creation.-
export_to_materials
(burnup_index, nuc_with_data=None) → openmc.material.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_CROSS_SECTIONS will be used.
Returns: mat_file – A modified Materials instance containing depleted material data and original isotopic compositions of non-depletable materials
Return type:
-
classmethod
from_hdf5
(filename)[source]¶ Load in depletion results from a previous file
Parameters: filename (str) – Path to depletion result file Returns: new – New instance of depletion results Return type: ResultsList
-
get_atoms
(mat, nuc, nuc_units='atoms', time_units='s')[source]¶ Get number of nuclides over time from a single material
Note
Initial values for some isotopes that do not appear in initial concentrations may be non-zero, depending on the value of
openmc.deplete.Operator
dilute_initial
. Theopenmc.deplete.Operator
adds isotopes according to this setting, which can be set to zero.Parameters: - mat (str) – Material name 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"}, optional) –
Units for the returned time array. Default is
"s"
to return the value in seconds.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_depletion_time
()[source]¶ Return an array of the average time to deplete a material
Note
Will have one fewer row than number of other methods, like
get_eigenvalues()
, because no depletion is performed at the final transport stageReturns: times – Vector of average time to deplete a single material across all processes and materials. Return type: numpy.ndarray
-
get_eigenvalue
()[source]¶ Evaluates the eigenvalue from a results list.
Returns: - times (numpy.ndarray) – Array of times in [s]
- eigenvalues (numpy.ndarray) – k-eigenvalue at each time. Column 0 contains the eigenvalue, while column 1 contains the associated uncertainty
-
get_reaction_rate
(mat, nuc, rx)[source]¶ Get reaction rate in a single material/nuclide over time
Note
Initial values for some isotopes that do not appear in initial concentrations may be non-zero, depending on the value of
openmc.deplete.Operator
dilute_initial
Theopenmc.deplete.Operator
adds isotopes according to this setting, which can be set to zero.Parameters: Returns: - times (numpy.ndarray) – Array of times in [s]
- rates (numpy.ndarray) – Array of reaction rates
-
get_step_where
(time, time_units='d', atol=1e-06, rtol=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
andrtol=math.inf
will return the closest index to the requested point.New in version 0.12.1.
Parameters: Returns: Return type:
-
get_times
(time_units='d') → numpy.ndarray[source]¶ Return the points in time that define the depletion schedule
New in version 0.12.1.
Parameters: time_units ({"s", "d", "h", "min"}, optional) – Return the vector in these units. Default is to convert to days Returns: 1-D vector of time points Return type: numpy.ndarray
-