openmc.deplete.SICELIIntegrator¶
- class openmc.deplete.SICELIIntegrator(operator: TransportOperator, timesteps: Sequence[float], power: Optional[Union[float, Sequence[float]]] = None, power_density: Optional[Union[float, Sequence[float]]] = None, source_rates: Optional[Sequence[float]] = None, timestep_units: str = 's', n_steps: int = 10, solver: str = 'cram48')[source]¶
Deplete using the SI-CE/LI CFQ4 algorithm.
Implements the stochastic implicit CE/LI predictor-corrector algorithm using the fourth order commutator-free integrator.
Detailed algorithm can be found in section 3.2 in Colin Josey’s thesis.
- Parameters
operator (openmc.deplete.abc.TransportOperator) – Operator to perform transport simulations
timesteps (iterable of float or iterable of tuple) – Array of timesteps. Note that values are not cumulative. The units are specified by the timestep_units argument when timesteps is an iterable of float. Alternatively, units can be specified for each step by passing an iterable of (value, unit) tuples.
power (float or iterable of float, optional) – Power of the reactor in [W]. A single value indicates that the power is constant over all timesteps. An iterable indicates potentially different power levels for each timestep. For a 2D problem, the power can be given in [W/cm] as long as the “volume” assigned to a depletion material is actually an area in [cm^2]. Either
power
,power_density
, orsource_rates
must be specified.power_density (float or iterable of float, optional) – Power density of the reactor in [W/gHM]. It is multiplied by initial heavy metal inventory to get total power if
power
is not specified.source_rates (float or iterable of float, optional) –
Source rate in [neutron/sec] or neutron flux in [neutron/s-cm^2] for each interval in
timesteps
New in version 0.12.1.
timestep_units ({'s', 'min', 'h', 'd', 'MWd/kg'}) – Units for values specified in the timesteps argument. ‘s’ means seconds, ‘min’ means minutes, ‘h’ means hours, and ‘MWd/kg’ indicates that the values are given in burnup (MW-d of energy deposited per kilogram of initial heavy metal).
n_steps (int, optional) – Number of stochastic iterations per depletion interval. Must be greater than zero. Default : 10
solver (str or callable, optional) –
If a string, must be the name of the solver responsible for solving the Bateman equations. Current options are:
cram16
- 16th order IPF CRAMcram48
- 48th order IPF CRAM [default]
If a function or other callable, must adhere to the requirements in
solver
.New in version 0.12.
- Variables
operator (openmc.deplete.abc.TransportOperator) – Operator to perform transport simulations
chain (openmc.deplete.Chain) – Depletion chain
timesteps (iterable of float) – Size of each depletion interval in [s]
power (iterable of float) – Power of the reactor in [W] for each interval in
timesteps
n_steps (int) – Number of stochastic iterations per depletion interval
solver (callable) –
Function that will solve the Bateman equations \(\frac{\partial}{\partial t}\vec{n} = A_i\vec{n}_i\) with a step size \(t_i\). Can be configured using the
solver
argument. User-supplied functions are expected to have the following signature:solver(A, n0, t) -> n1
whereA
is ascipy.sparse.csc_matrix
making up the depletion matrixn0
is a 1-Dnumpy.ndarray
of initial compositions for a given material in atoms/cm3t
is a float of the time step size in seconds, andn1
is anumpy.ndarray
of compositions at the next time step. Expected to be of the same shape asn0
New in version 0.12.
- __call__(n_bos, bos_rates, dt, source_rate, _i=None)[source]¶
Perform the integration across one time step
- Parameters
n_bos (list of numpy.ndarray) – List of atom number arrays for each material. Each array in the list contains the number of [atom] of each nuclide.
bos_rates (openmc.deplete.ReactionRates) – Reaction rates from operator
dt (float) – Time in [s] for the entire depletion interval
source_rate (float) – Power in [W] or source rate in [neutron/sec]
_i (int, optional) – Current depletion step index. Not used
- Returns
proc_time (float) – Time spent in CRAM routines for all materials in [s]
n_bos_list (list of list of numpy.ndarray) – Concentrations at each of the intermediate points with the final concentration as the last element
op_results (list of openmc.deplete.OperatorResult) – Eigenvalue and reaction rates from intermediate transport simulations
- __iter__()¶
Return pair of time step in [s] and source rate in [W] or [neutron/sec]
- __len__()¶
Return integer number of depletion intervals
- add_transfer_rate(material: Union[str, int, Material], components: Sequence[str], transfer_rate: float, transfer_rate_units: str = '1/s', destination_material: Optional[Union[str, int, Material]] = None)¶
Add transfer rates to depletable material.
- Parameters
material (openmc.Material or str or int) – Depletable material
components (list of str) – List of strings of elements and/or nuclides that share transfer rate. A transfer rate for a nuclide cannot be added to a material alongside a transfer rate for its element and vice versa.
transfer_rate (float) – Rate at which elements are transferred. A positive or negative values set removal of feed rates, respectively.
destination_material (openmc.Material or str or int, Optional) – Destination material to where nuclides get fed.
transfer_rate_units ({'1/s', '1/min', '1/h', '1/d', '1/a'}) – Units for values specified in the transfer_rate argument. ‘s’ means seconds, ‘min’ means minutes, ‘h’ means hours, ‘a’ means Julian years.
- integrate(output: bool = True, path: Union[str, PathLike] = 'depletion_results.h5')¶
Perform the entire depletion process across all steps
- Parameters
output (bool, optional) – Indicate whether to display information about progress
path (PathLike) –
Path to file to write. Defaults to ‘depletion_results.h5’.
New in version 0.15.0.