openmc.deplete.CELIIntegrator¶
-
class
openmc.deplete.
CELIIntegrator
(operator, timesteps, power=None, power_density=None, timestep_units='s', solver='cram48')[source]¶ Deplete using the CE/LI CFQ4 algorithm.
Implements the CE/LI Predictor-Corrector algorithm using the fourth order commutator-free integrator.
“CE/LI” stands for constant extrapolation on predictor and linear interpolation on corrector. This algorithm is mathematically defined as:
\[\begin{aligned} y' &= A(y, t) y(t) \\ A_0 &= A(y_n, t_n) \\ y_p &= \text{expm}(h A_0) y_n \\ A_1 &= A(y_p, t_n + h) \\ y_{n+1} &= \text{expm}(\frac{h}{12} A_0 + \frac{5h}{12} A1) \text{expm}(\frac{5h}{12} A_0 + \frac{h}{12} A1) y_n \end{aligned} \]Parameters: - operator (openmc.deplete.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
orpower_density
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 speficied. - 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).
New in version 0.12.
- 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.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
- 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.csr_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__
(bos_conc, rates, dt, power, _i=None)[source]¶ Perform the integration across one time step
Parameters: - bos_conc (numpy.ndarray) – Initial concentrations for all nuclides in [atom]
- rates (openmc.deplete.ReactionRates) – Reaction rates from operator
- dt (float) – Time in [s] for the entire depletion interval
- power (float) – Power of the system in [W]
- _i (int, optional) – Current iteration count. Not used
Returns: - proc_time (float) – Time spent in CRAM routines for all materials in [s]
- conc_list (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 simulation
-
__iter__
()¶ Return pairs of time steps in [s] and powers in [W]
-
__len__
()¶ Return integer number of depletion intervals
-
integrate
()¶ Perform the entire depletion process across all steps