9. Iterative Methods
9.1. Overview
OpenSn’s linear Boltzmann solver uses several levels of iteration. The most important distinction for users is that not all iteration controls live in the same place.
There are three levels to keep in mind:
inner groupset solves
across-groupset iterations
solver-specific outer iterations, such as k-eigenvalue iterations
In an OpenSn input, the inner groupset solver is selected in each groupset block,
while across-groupset controls are set in the problem options block.
Example:
phys = DiscreteOrdinatesProblem(
mesh=grid,
num_groups=20,
groupsets=[
{
"groups_from_to": (0, 9),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 200,
"gmres_restart_interval": 30,
},
{
"groups_from_to": (10, 19),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 200,
"gmres_restart_interval": 30,
},
],
xs_map=[{"block_ids": [0], "xs": xs}],
options={
"max_ags_iterations": 100,
"ags_tolerance": 1.0e-6,
"ags_convergence_check": "l2",
},
)
For most users, the practical questions are:
which inner solver should be used in each groupset?
should DSA be enabled?
does the problem need multiple groupsets and AGS control?
is the top-level solver a steady-state solve, a transient solve, a power iteration k-eigen solve, or a nonlinear k-eigen solve?
9.2. Where the Iteration Controls Live
9.2.1. Groupset controls
These are specified inside each entry of groupsets=[...]:
inner_linear_methodl_abs_toll_max_itsgmres_restart_intervalallow_cyclesapply_wgdsawgdsa_l_abs_tolwgdsa_l_max_itswgdsa_verbosewgdsa_petsc_optionsapply_tgdsatgdsa_l_abs_toltgdsa_l_max_itstgdsa_verbosetgdsa_petsc_options
These settings control the solve performed within a single groupset.
9.2.2. Problem-level controls
These are specified in the problem options block:
max_ags_iterationsags_toleranceags_convergence_checkverbose_inner_iterationsverbose_ags_iterationsverbose_outer_iterations
These settings control how the solver coordinates multiple groupsets and how much iteration information is printed.
9.2.3. Solver-level controls
Some solvers add their own iteration settings.
For pyopensn.solver.PowerIterationKEigenSolver:
max_itersk_tolreset_phi0
For pyopensn.solver.NonLinearKEigenSolver:
nl_abs_tolnl_rel_tolnl_sol_tolnl_max_itsl_abs_toll_rel_toll_div_toll_max_itsl_gmres_restart_intvll_gmres_breakdown_tolreset_phi0num_initial_power_iterations
9.3. Inner Groupset Methods
The available groupset inner methods are:
"classic_richardson""petsc_richardson""petsc_gmres""petsc_bicgstab"
These are selected with inner_linear_method in each groupset.
Example:
groupset = {
"groups_from_to": (0, 9),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 300,
"gmres_restart_interval": 50,
}
9.3.1. What these methods are solving
For a discrete ordinates groupset, the inner iteration repeatedly applies the inverse transport sweep operator and updates the groupset flux moments until the chosen convergence criterion is met.
The main difference between the available methods is how that fixed-point or linear solve is managed:
classic_richardsonperforms an explicit source iteration loop in OpenSn codethe
petsc_*methods use PETSc Krylov solvers wrapped around the sweep operator
As a practical rule, PETSc methods are usually the better starting point for
difficult problems. classic_richardson is simple and useful, but it is
generally the least robust option on strongly scattering or otherwise difficult
systems unless acceleration is also enabled.
9.3.2. classic_richardson
This is OpenSn’s explicit source-iteration style method. On each iteration it:
rebuilds the groupset source from the current iterate
applies the inverse transport sweep
optionally applies WGDSA and TGDSA corrections
checks convergence using pointwise changes in
phiand delayed angular flux
Why use it:
simple algorithmic behavior
useful for debugging and method studies
pairs naturally with DSA for source-iteration style solves
by far the least memory-intensive of the groupset inner solver methods
Why not use it by default:
usually slower and less robust than Krylov methods on harder problems
convergence can deteriorate badly for optically thick or highly scattering cases without acceleration
Example:
{
"groups_from_to": (0, 0),
"angular_quadrature": quad,
"inner_linear_method": "classic_richardson",
"l_abs_tol": 1.0e-8,
"l_max_its": 200,
}
9.3.3. petsc_richardson
This uses PETSc’s Richardson iteration around the transport operator.
Why use it:
useful when a simple stationary iteration is desired but with PETSc-managed setup and convergence handling
A special case worth knowing:
if ``inner_linear_method=”petsc_richardson”`` and ``l_max_its=1``, the code effectively performs a single sweep-style update
For most problems, this is not the first PETSc method to try. GMRES is usually the stronger default.
9.3.4. petsc_gmres
This uses PETSc GMRES around the transport operator.
For most steady-state transport problems, this is the best general-purpose starting choice.
Why use it:
usually the most robust of the groupset solver methods
often the best choice for difficult scattering-dominated problems
works with DSA preconditioning
Relevant control:
gmres_restart_intervalcontrols how many Krylov vectors GMRES keeps before it restarts
Warning
GMRES can incur substantial memory usage because it stores restart vectors.
Increasing gmres_restart_interval increases the Krylov storage footprint,
so use larger values only when the convergence benefit justifies the extra
memory.
Example:
{
"groups_from_to": (0, 19),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 300,
"gmres_restart_interval": 30,
}
9.3.5. petsc_bicgstab
This uses PETSc BiCGStab around the transport operator.
Why use it:
as an alternative to GMRES when memory or restart behavior is undesirable
as a solver-choice experiment for problems where GMRES is not ideal
In most common workflows, GMRES is still the more typical choice.
9.4. Convergence Controls
9.4.1. Groupset tolerance and iteration count
Each groupset has:
l_abs_tol: default1.0e-6l_max_its: default200
These control the inner solve stopping behavior.
For PETSc-based groupset methods, the residual printed in the log is the scaled residual used by OpenSn’s custom convergence test for WGS solves.
For classic_richardson, the convergence test is based on pointwise changes
in the solution, together with a spectral-radius estimate.
9.4.1.1. l_abs_tol
This is the main convergence tolerance for the groupset inner solve.
What it tests:
for PETSc-based groupset methods, it tests the scaled WGS residual used by the OpenSn convergence check
for
classic_richardson, it is the stopping tolerance for the pointwise flux-change test used by the Richardson iteration machinery
Reasonable values:
1.0e-6is a good default for most production runs1.0e-4to1.0e-5is often acceptable for setup, debugging, or quick parameter studies1.0e-7to1.0e-8is reasonable when eigenvalues, localized responses, or balance quantities need tighter transport convergence
If the solve repeatedly reaches the outer iterations while the inner residual is still large, the inner tolerance is probably too loose for that problem.
9.4.1.2. l_max_its
This is the hard cap on the number of groupset inner iterations.
What it tests:
nothing by itself; it is a safeguard that stops the inner solve if the chosen convergence test is not met in time
Reasonable values:
200is a good default100to300is a normal range for most problemshigher values can be justified for very difficult scattering-dominated problems, but routinely hitting the cap usually means the method, groupset structure, or acceleration setup needs attention
9.4.2. What tolerance to choose
As a practical rule:
use looser tolerances while setting up and debugging a model
tighten tolerances when k-eigenvalues, localized responses, or balance quantities need to be resolved accurately
if outer iterations are still moving substantially, over-solving the inner iterations may not be the best use of work
This is especially important in multi-groupset and eigenvalue problems, where the inner solve is only one part of the total iteration process.
9.4.3. GMRES restart interval
gmres_restart_interval only matters for petsc_gmres.
Smaller restart values use less Krylov storage but may slow or degrade convergence. Larger restart values can improve robustness but cost more memory and work per restart cycle.
What it tests:
it does not change the convergence test directly
it controls how many GMRES Krylov vectors are retained before restart
Reasonable values:
30is a good default20to50is a reasonable first rangelarger values such as
60to100can help difficult problems if the extra memory is acceptable
Note
On large production cases, the memory cost of GMRES restart vectors can be a
practical limiter. If memory grows too much, reduce
gmres_restart_interval first before assuming the only fix is more
hardware.
9.4.4. allow_cycles
allow_cycles controls whether sweep cycles are allowed in the groupset.
For most users, this is not the first iterative-method parameter to tune, but it is part of the groupset solve setup and can matter on meshes whose sweep dependencies contain cycles.
9.5. Across-Groupset Iteration
When a problem has multiple groupsets, OpenSn may perform across-groupset iterations, usually abbreviated AGS.
The relevant problem options are:
max_ags_iterationsags_toleranceags_convergence_checkverbose_ags_iterations
9.5.1. What AGS does
AGS coordinates repeated solves of the groupsets until the full multigroup solution is consistent across groupset boundaries.
In OpenSn, AGS:
executes each WGS solver in sequence
compares the new and previous full
phiiterateuses either an
l2orpointwiseconvergence checkstops when the AGS tolerance is satisfied or the iteration limit is reached
Why this matters:
if the problem has only one groupset, AGS is irrelevant
if the problem has multiple groupsets, AGS is part of the main solve algorithm, not just a logging detail
9.5.2. ags_convergence_check
Allowed values are:
"l2""pointwise"
l2 is the default and is usually the simplest first choice.
pointwise is stricter in a different sense: it tracks the maximum local
relative change in the flux vector. It can be useful when localized solution
changes matter more than the global norm.
9.5.3. ags_tolerance
This is the stopping tolerance for the AGS iteration.
What it tests:
the change between successive full multigroup
phiiterates across AGS passes, using either thel2orpointwisecheck selected byags_convergence_check
Reasonable values:
1.0e-6is a good default1.0e-5is often sufficient for debugging or coarse studies1.0e-7to1.0e-8may be appropriate when strong cross-groupset coupling must be converged tightly
In general, the groupset inner solves should be at least as tight as the AGS solve, and usually somewhat tighter.
9.5.4. max_ags_iterations
This is the hard cap on AGS passes.
What it tests:
nothing by itself; it simply limits the total number of across-groupset iterations
Reasonable values:
100is a good default20to50is often enough for easier multi-groupset problemsif this limit is reached regularly, the real issue is usually strong cross-groupset coupling, an overly aggressive groupset split, or the need for tighter inner solves
Example:
options = {
"max_ags_iterations": 100,
"ags_tolerance": 1.0e-6,
"ags_convergence_check": "l2",
"verbose_ags_iterations": True,
}
9.6. DSA Options
The groupset-level acceleration toggles are:
apply_wgdsaapply_tgdsa
with associated controls:
wgdsa_l_abs_tolwgdsa_l_max_itswgdsa_verbosewgdsa_petsc_optionstgdsa_l_abs_toltgdsa_l_max_itstgdsa_verbosetgdsa_petsc_options
9.6.1. WGDSA
WGDSA stands for within-group diffusion synthetic acceleration. Enabling
apply_wgdsa creates a diffusion MIP solve over the current groupset and
applies the resulting correction to the updated flux.
Why use it:
to accelerate scattering-dominated inner iterations
especially useful with Richardson-style iteration
can also help PETSc-based inner solves when used as part of the transport solve workflow
WGDSA controls:
wgdsa_l_abs_tol: stopping tolerance for the diffusion correction solvewgdsa_l_max_its: hard cap on WGDSA iterations
Reasonable values:
start with
wgdsa_l_abs_tol=1.0e-4andwgdsa_l_max_its=30if WGDSA is clearly helping but not converging tightly enough, reduce the tolerance to
1.0e-5or raise the iteration cap modestlyit is usually not necessary to solve the DSA system more tightly than the transport problem itself
9.6.2. TGDSA
TGDSA stands for two-grid diffusion synthetic acceleration. Enabling
apply_tgdsa constructs a collapsed one-group diffusion correction using
two-grid information derived from the groupset cross sections, then projects the
correction back to the groupset spectrum.
TGDSA controls:
tgdsa_l_abs_tol: stopping tolerance for the two-grid diffusion solvetgdsa_l_max_its: hard cap on TGDSA iterations
Reasonable values:
the same starting values as WGDSA are usually sensible:
tgdsa_l_abs_tol=1.0e-4andtgdsa_l_max_its=30TGDSA is a correction solve, so very tight tolerances are usually unnecessary
Why use it:
to accelerate multigroup error modes that are not well treated by within-group acceleration alone
9.6.3. Practical guidance
For many problems:
start without DSA
if convergence is too slow, try WGDSA first
add TGDSA when multigroup coupling remains difficult
DSA is most useful when there is an actual iteration problem to solve. It is not automatically beneficial for every easy case.
Example:
{
"groups_from_to": (0, 39),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 200,
"apply_wgdsa": True,
"wgdsa_l_abs_tol": 1.0e-4,
"wgdsa_l_max_its": 30,
"wgdsa_verbose": False,
}
9.7. Solver-Specific Outer Iteration
9.7.1. Steady-state source solves
For a standard steady-state source problem, the main iteration controls are the groupset inner method and, when multiple groupsets are present, AGS.
Example:
ss_solver = SteadyStateSourceSolver(problem=phys)
ss_solver.Initialize()
ss_solver.Execute()
9.7.2. Transient solves
For transient problems, each timestep calls into the same underlying transport iteration machinery. That means the groupset inner method and AGS settings still matter, but now they matter per timestep.
If a transient run is unexpectedly expensive, the first iterative-method settings to examine are usually:
the chosen groupset inner method
the inner tolerance
whether angular flux storage is enabled for the transient formulation
9.7.3. Power iteration k-eigen solve
pyopensn.solver.PowerIterationKEigenSolver adds its own outer
iteration on top of the transport solve.
Relevant parameters:
max_itersk_tolreset_phi0
What it does:
sets the fission source from the previous iterate
performs the transport solve through AGS and the groupset inner solves
updates
k_effrepeats until the eigenvalue change is below
k_tolormax_itersis reached
k_tol is the stopping tolerance on the eigenvalue change between successive
power iterations. 1.0e-8 to 1.0e-10 is a reasonable production range,
with 1.0e-10 being common when a fairly tight k_eff is desired.
max_iters is the hard cap on power iterations. 500 to 1000 is a
reasonable starting range. If that cap is reached often, the issue is usually
not the cap itself but the overall transport convergence strategy.
Example:
keigen = PowerIterationKEigenSolver(
problem=phys,
max_iters=500,
k_tol=1.0e-10,
)
keigen.Initialize()
keigen.Execute()
9.7.4. Nonlinear k-eigen solve
pyopensn.solver.NonLinearKEigenSolver uses a separate nonlinear
solve path.
Relevant parameters:
nonlinear tolerances:
nl_abs_tol,nl_rel_tol,nl_sol_tol,nl_max_itslinear tolerances inside the nonlinear solve:
l_abs_tol,l_rel_tol,l_div_tol,l_max_itsGMRES controls for the nonlinear linearization solves:
l_gmres_restart_intvlandl_gmres_breakdown_toloptional warm start via
num_initial_power_iterations
Why use it:
when the nonlinear eigenvalue solve is preferred over plain power iteration
when convergence behavior of the nonlinear method is better suited to the problem
This is a different top-level algorithm from groupset WGS iteration. Its
internal linear solves are controlled by the l_* parameters on
NonLinearKEigenSolver, not by the groupset inner_linear_method or the
groupset l_* settings.
Important
This is different from the other solver paths discussed in this section.
NonLinearKEigenSolver uses its own internal linear solve. If a nonlinear
k-eigen run needs tuning, start with the solver-level nl_* and l_*
parameters on pyopensn.solver.NonLinearKEigenSolver, not with
the groupset inner-solver settings.
What the main nonlinear controls do:
nl_abs_tol: absolute nonlinear residual tolerancenl_rel_tol: relative nonlinear residual tolerancenl_sol_tol: solution-update tolerance for the nonlinear solvenl_max_its: hard cap on nonlinear iterations
Reasonable values:
1.0e-8is a sensible starting point fornl_abs_tolandnl_rel_tol25to50is a reasonable starting range fornl_max_itsnl_sol_tolis usually left at its default unless there is a specific reason to tune the solution-update criterion
What the internal linear controls do:
l_abs_tol: absolute tolerance for the linear solve inside each nonlinear stepl_rel_tol: relative tolerance for that same linear solvel_div_tol: divergence guard for the linear solvel_max_its: hard cap on linear iterations inside each nonlinear stepl_gmres_restart_intvl: restart size for the internal GMRES solvel_gmres_breakdown_tol: breakdown guard for that GMRES solve
Reasonable values:
1.0e-8is a reasonable starting value forl_abs_tolandl_rel_tol50is a reasonable first value forl_max_its30is a good default forl_gmres_restart_intvll_div_tolandl_gmres_breakdown_tolare usually left at their defaults unless the nonlinear linearization is clearly running into solver pathologies
9.8. Choosing a Method
For most users, a sensible starting progression is:
Start with
petsc_gmres.Use the default AGS settings unless there is a clear need to change them.
If convergence is slow on scattering-dominated problems, enable WGDSA.
If multigroup convergence remains difficult, consider TGDSA and/or revisit groupset splitting.
Use
classic_richardsonwhen a simple source-iteration style method is specifically desired.
Keep in mind that petsc_gmres is not memory-free: the restart interval sets
how many Krylov vectors are retained between restarts, so aggressive values can
materially increase memory usage.
If the problem is easy, these choices may not matter much. If the problem is difficult, the best improvements usually come from:
a better inner method
appropriate DSA
sensible groupset decomposition
tolerances that match the real accuracy target
9.9. Examples
9.9.1. Single groupset with GMRES
groupsets = [
{
"groups_from_to": (0, 9),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 200,
"gmres_restart_interval": 30,
}
]
9.9.2. Classic Richardson with WGDSA
groupsets = [
{
"groups_from_to": (0, 19),
"angular_quadrature": quad,
"inner_linear_method": "classic_richardson",
"l_abs_tol": 1.0e-8,
"l_max_its": 200,
"apply_wgdsa": True,
"wgdsa_l_abs_tol": 1.0e-4,
"wgdsa_l_max_its": 30,
}
]
9.9.3. Two groupsets with AGS control
phys = DiscreteOrdinatesProblem(
mesh=grid,
num_groups=40,
groupsets=[
{
"groups_from_to": (0, 19),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 200,
},
{
"groups_from_to": (20, 39),
"angular_quadrature": quad,
"inner_linear_method": "petsc_gmres",
"l_abs_tol": 1.0e-6,
"l_max_its": 200,
},
],
xs_map=[{"block_ids": [0], "xs": xs}],
options={
"max_ags_iterations": 100,
"ags_tolerance": 1.0e-6,
"ags_convergence_check": "pointwise",
"verbose_ags_iterations": True,
},
)
9.9.4. Power iteration k-eigen solve
keigen = PowerIterationKEigenSolver(
problem=phys,
max_iters=1000,
k_tol=1.0e-10,
reset_phi0=True,
)
9.9.5. Nonlinear k-eigen solve
keigen = NonLinearKEigenSolver(
problem=phys,
nl_abs_tol=1.0e-8,
nl_rel_tol=1.0e-8,
nl_max_its=50,
l_abs_tol=1.0e-8,
l_rel_tol=1.0e-8,
l_max_its=50,
l_gmres_restart_intvl=30,
num_initial_power_iterations=5,
)
9.10. Recommendations
Use
petsc_gmresas the default first choice unless there is a specific reason to prefer another method.Remember that GMRES stores restart vectors, so aggressive restart sizes can significantly increase memory usage.
Use
classic_richardsonwhen a simple source-iteration style solve is desired, especially for studies and debugging.Treat AGS settings as important only when the problem has multiple groupsets.
Try WGDSA before TGDSA when the main issue is slow scattering convergence.
Tighten tolerances only as far as the quantities of interest actually require.
Revisit the iterative-method settings when changing problem difficulty, groupset structure, or solver type. A choice that is good for a fixed-source problem may not be the best one for a transient or k-eigen solve.