openmc.deplete.LEQIIntegrator

class openmc.deplete.LEQIIntegrator(operator, timesteps, power=None, power_density=None, source_rates=None, timestep_units='s', solver='cram48')[source]

Deplete using the LE/QI CFQ4 algorithm.

Implements the LE/QI Predictor-Corrector algorithm using the fourth order commutator-free integrator.

“LE/QI” stands for linear extrapolation on predictor and quadratic interpolation on corrector. This algorithm is mathematically defined as:

\[\begin{aligned} \mathbf{A}_{-1} &= \mathbf{A}(\mathbf{n}_{i-1}) \\ \mathbf{A}_0 &= \mathbf{A}(\mathbf{n}_i) \\ \mathbf{F}_1 &= \frac{-h_i}{12h_{i-1}} \mathbf{A}_{-1} + \frac{6h_{i-1} + h_i}{12h_{i-1}} \mathbf{A}_0 \\ \mathbf{F}_2 &= \frac{-5h_i}{12h_{i-1}} \mathbf{A}_{-1} + \frac{6h_{i-1} + 5h_i}{12h_{i-1}} \mathbf{A}_0 \\ \mathbf{n}_{i+1}^p &= \exp (h_i \mathbf{F}_1) \exp(h_i \mathbf{F}_2) \mathbf{n}_i \\ \mathbf{A}_1 &= \mathbf{A}(\mathbf{n}_{i+1}^p) \\ \mathbf{F}_3 &= \frac{-h_i^2}{12 h_{i-1} (h_{i-1} + h_i)} \mathbf{A}_{-1} + \frac{5 h_{i-1}^2 + 6 h_i h_{i-1} + h_i^2}{12 h_{i-1} (h_{i-1} + h_i)} \mathbf{A}_0 + \frac{h_{i-1}}{12 (h_{i-1} + h_i)} \mathbf{A}_1 \\ \mathbf{F}_4 &= \frac{-h_i^2}{12 h_{i-1} (h_{i-1} + h_i)} \mathbf{A}_{-1} + \frac{h_{i-1}^2 + 2 h_i h_{i-1} + h_i^2}{12 h_{i-1} (h_{i-1} + h_i)} \mathbf{A}_0 + \frac{5 h_{i-1}^2 + 4 h_i h_{i-1}}{12 h_{i-1} (h_{i-1} + h_i)} \mathbf{A}_1 \\ \mathbf{n}_{i+1} &= \exp(h_i \mathbf{F}_4) \exp(h_i \mathbf{F}_3) \mathbf{n}_i \end{aligned} \]

It is initialized using the CE/LI algorithm.

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, or source_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', 'a', 'MWd/kg'}) – Units for values specified in the timesteps argument. ‘s’ means seconds, ‘min’ means minutes, ‘h’ means hours, ‘a’ means Julian years and ‘MWd/kg’ indicates that the values are given in burnup (MW-d of energy deposited per kilogram of initial heavy metal).

  • 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 CRAM

    • cram48 - 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]

  • source_rates (iterable of float) – Source rate in [W] or [neutron/sec] 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 where

    • A is a scipy.sparse.csr_matrix making up the depletion matrix

    • n0 is a 1-D numpy.ndarray of initial compositions for a given material in atoms/cm3

    • t is a float of the time step size in seconds, and

    • n1 is a numpy.ndarray of compositions at the next time step. Expected to be of the same shape as n0

    New in version 0.12.

__call__(bos_conc, bos_rates, dt, source_rate, i)[source]

Perform the integration across one time step

Parameters
  • 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

  • source_rate (float) – Power in [W] or source rate in [neutron/sec]

  • i (int) – Current depletion step index

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 pair of time step in [s] and source rate in [W] or [neutron/sec]

__len__()

Return integer number of depletion intervals

integrate(final_step=True, output=True)

Perform the entire depletion process across all steps

Parameters
  • final_step (bool, optional) –

    Indicate whether or not a transport solve should be run at the end of the last timestep.

    New in version 0.12.1.

  • output (bool, optional) –

    Indicate whether to display information about progress

    New in version 0.13.1.