openmc.deplete.CF4Integrator

class openmc.deplete.CF4Integrator(operator, timesteps, power=None, power_density=None)[source]

Deplete using the CF4 algorithm.

Implements the fourth order commutator-free Lie algorithm. This algorithm is mathematically defined as:

\[\begin{aligned} F_1 &= h A(y_0) \\ y_1 &= \text{expm}(1/2 F_1) y_0 \\ F_2 &= h A(y_1) \\ y_2 &= \text{expm}(1/2 F_2) y_0 \\ F_3 &= h A(y_2) \\ y_3 &= \text{expm}(-1/2 F_1 + F_3) y_1 \\ F_4 &= h A(y_3) \\ y_4 &= \text{expm}( 1/4 F_1 + 1/6 F_2 + 1/6 F_3 - 1/12 F_4) \text{expm}(-1/12 F_1 + 1/6 F_2 + 1/6 F_3 + 1/4 F_4) y_0 \end{aligned} \]
Parameters:
  • operator (openmc.deplete.TransportOperator) – Operator to perform transport simulations
  • timesteps (iterable of float) – Array of timesteps in units of [s]. Note that values are not cumulative.
  • 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 or power_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.
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
__call__(bos_conc, bos_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]
  • bos_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 depletion step index. 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 simulations

__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