pyopensn.solver.NonLinearKEigenSolver

class pyopensn.solver.NonLinearKEigenSolver

Non-linear k-eigenvalue solver.

Wrapper of opensn::NonLinearKEigenSolver.

Supports one or more groupsets when the groupsets run without WGDSA/TGDSA. If WGDSA or TGDSA is enabled on any groupset, the nonlinear k-eigen solve must use a single groupset.

Advance(self: pyopensn.solver.Solver) None

Advance time values function.

ComputeBalanceTable(self: pyopensn.solver.NonLinearKEigenSolver) dict

Compute and return the global balance table using the solver’s normalization. This is a collective operation and must be called on all ranks.

Returns:

Dictionary with the following entries:

  • absorption_rate: Global absorption rate, approximately integral sigma_a * phi dV summed over groups and the full domain.

  • production_rate: Global volumetric production/source rate used by the solver, approximately integral Q dV summed over groups and the full domain.

  • inflow_rate: Global incoming boundary contribution integrated over incoming angular flux on boundaries.

  • outflow_rate: Global outgoing boundary contribution accumulated from face outflow tallies.

  • balance: Rate balance, production_rate + inflow_rate - absorption_rate - outflow_rate.

Return type:

dict

Notes

For k-eigenvalue balance reporting, this solver scales the production term by 1 / k_eff before forming both production_rate and balance.

Execute(self: pyopensn.solver.Solver) None

Execute the solver.

GetEigenvalue(self: pyopensn.solver.NonLinearKEigenSolver) float

Return the current k-eigenvalue.

Initialize(self: pyopensn.solver.Solver) None

Initialize the solver.

__init__(self: pyopensn.solver.NonLinearKEigenSolver, **kwargs) None

Construct a non-linear k-eigenvalue solver.

Parameters:
  • problem (pyopensn.solver.DiscreteOrdinatesProblem) – Existing discrete ordinates problem instance. Multiple groupsets are supported when groupset WGDSA/TGDSA is disabled. If WGDSA or TGDSA is enabled on any groupset, the nonlinear k-eigen solve must use a single groupset.

  • nl_abs_tol (float, default=1.0e-8) – Non-linear absolute tolerance.

  • nl_rel_tol (float, default=1.0e-8) – Non-linear relative tolerance.

  • nl_sol_tol (float, default=1.0e-50) – Non-linear solution tolerance.

  • nl_max_its (int, default=50) – Non-linear algorithm maximum iterations.

  • l_abs_tol (float, default=1.0e-8) – Linear absolute tolerance.

  • l_rel_tol (float, default=1.0e-8) – Linear relative tolerance.

  • l_div_tol (float, default=1.0e6) – Linear divergence tolerance.

  • l_max_its (int, default=50) – Linear algorithm maximum iterations.

  • l_gmres_restart_intvl (int, default=30) – GMRES restart interval.

  • l_gmres_breakdown_tol (float, default=0.1) – GMRES breakdown tolerance.

  • reset_phi0 (bool, default=True) – If true, reinitializes scalar fluxes to 1.0.

  • num_initial_power_iterations (int, default=0) – Number of initial power iterations before the non-linear solve.

Notes

PETSc convergence failures and iteration limits are reported through the solver status and log output, consistent with the other transport solvers. Invalid residual states are reported to PETSc as function-domain errors so SNES can backtrack or terminate with the appropriate reason.