openmc.deplete.cram.IPFCramSolver
- class openmc.deplete.cram.IPFCramSolver(alpha, theta, alpha0)[source]
CRAM depletion solver that uses incomplete partial factorization
Provides a
__call__()that utilizes an incomplete partial factorization (IPF) for the Chebyshev Rational Approximation Method (CRAM), as described in the following paper: M. Pusa, “Higher-Order Chebyshev Rational Approximation Method and Application to Burnup Equations,” Nucl. Sci. Eng., 182:3, 297-318.When substeps > 1, the time interval is split into substeps identical sub-intervals and LU factorizations are reused across them, as described in: A. Isotalo and M. Pusa, “Improving the Accuracy of the Chebyshev Rational Approximation Method Using Substeps,” Nucl. Sci. Eng., 183:1, 65-77.
- Parameters:
alpha (numpy.ndarray) – Complex residues of poles used in the factorization. Must be a vector with even number of items.
theta (numpy.ndarray) – Complex poles. Must have an equal size as
alpha.alpha0 (float) – Limit of the approximation at infinity
- Variables:
alpha (numpy.ndarray) – Complex residues of poles
thetain the incomplete partial factorization. Denoted as \(\tilde{\alpha}\)theta (numpy.ndarray) – Complex poles \(\theta\) of the rational approximation
alpha0 (float) – Limit of the approximation at infinity
- __call__(A, n0, dt, substeps=1)[source]
Solve depletion equations using IPF CRAM
- Parameters:
A (scipy.sparse.csc_array) – Sparse transmutation matrix
A[j, i]desribing rates at which isotopeitransmutes to isotopejn0 (numpy.ndarray) – Initial compositions, typically given in number of atoms in some material or an atom density
dt (float) – Time [s] of the specific interval to be solved
substeps (int, optional) – Number of substeps per depletion interval.
- Returns:
Final compositions after
dt- Return type: