Steady

Executioner for steady-state simulations.

Overview

Steady is a general solver for discrete steady-state nonlinear or linear problem:

(1)

By default the line-search Newton in PETSc is used with the PJFNK (preconditioned Jacobian-free Newton Krylov) method. At each Newton iteration the executioner solves

(2)

where

is the Jacobian matrix evaluated at . Jacobian matrix depends on for general nonlinear problems while it is constant for linear problems. The right hand side is also typically referred to as the residual at of Eq. (1). The Krylov methods are employed for solving the above linear equation, which requires only the evaluation of the matrix-vector product . Within the MOOSE framework, the Jacobian-Free Newton Krylov method is used that approximates matrix vector products by the Finite-Difference like approximation:

(3)

where the scalar value is chosen by PETSc automatically to approximate accurately for the linear solve. It is noted that for a linear problem of which can be expressed as , where matrix is the Jacobian independent on and is the right-hand-side vector, the right hand side of Eq. (3) is independent on . Section 5.5 of PETSc user's manual on matrix-free methods details the algorithm for choosing the value of . It is actually the PETSc option -mat_mffd_err controls the but not -snes_mf_err unless we set -snes_mf_version to 2 other than the default 1. This could be changed in future PETSc updates.

The Krylov methods typically also require an approximation of the actual Jacobian for pre-conditioning the Krylov solution at each linear iteration. Note, the preconditioning matrix is seldom the exact Jacobian because it would require too much computational time and memory to compute, and in some cases is simply impossible to compute. By default the type of Krylov method in use is GMRES because it does not have assumptions on the underlying Jacobian. The initial guess for each linear solve is always set to zero, which implies that the initial linear residual is the same of the nonlinear residual. The residual norm at each linear iteration is evaluated by PETSc, for instance, during updating the Hessenberg matrix if GMRES method is used. At the conclusion of the nonlinear iteration, the solution is updated as follows

(4)

where is determined by the line-search algorithm. We can see that at each nonlinear or Newton iteration, we will need to update the preconditioning matrix and evaluate the residual with the updated solution. At each linear iteration, we simply need a residual evaluation and the operation of the preconditioner built from the preconditioning matrix. In PETSc the preconditioner type refers to the method to obtain an approximation of the inverse of and not a means to compute the elements of . It is noted that the default preconditioner type depends on the number of processors and also depends on the assembled preconditioning matrix . Typically incomplete LU (PCILU) is the default type with one processor and block Jacobi (PCBJACOBI) is the type with multiple processors. Consequently, you will not see the same convergence with the different number of processors. Note, there are two approximations in play here: (1) the Jacobian is approximated by a matrix that is easier to compute, and (2) the matrix is inverted approximately. The preconditioning matrix can be viewed with the PETSc option -ksp_view_pmat.

Solve Type

The general method in which the nonlinear system is solved is controlled by the "solve_type" parameter. Below is a description of each of the options:

  • PJFNK is the default solve type. It makes the executioner perform Jacobian-free linear solves at each Newton iteration with the preconditioner built from the preconditioning matrix . By default, the preconditioning matrix is block-diagonal with each block corresponding to a single MOOSE variable without custom preconditioning, refer to Preconditioning. Off-diagonal Jacobian terms are ignored. It essentially activates the matrix-free Jacobian-vector products, and the preconditioning matrix.

  • JFNK means there is no preconditioning during the Krylov solve. No Jacobian will be assembled. It essentially activates the matrix-free Jacobian-vector products and no preconditioning matrix.

  • LINEAR will use PETSc control parameter -snes_type ksponly to set the type of SNES for solving the linear system. Note that it only works when you have an exact Jacobian because it is not activating matrix-free calculations.

  • NEWTON means PETSc will use the Jacobian provided by kernels (typically not exact) to do the Krylov solve. If the Jacobian is not exact, Newton update in Eq. (4) will not reduce the residual effectively and typically results into an unconverged Newton iteration.

  • FD means the Jacobian is assembled via a finite differencing method. This is costly and should used only for testing purpose.

Preconditioning

  • Krylov methods need preconditioning to be efficient (or even effective!).

  • Even though the Jacobian is never formed, JFNK methods still require preconditioning.

  • MOOSE's automatic (without user intervention) preconditioning is fairly minimal.

  • Many options exist for implementing improved preconditioning in MOOSE.

Preconditioned JFNK

  • Using right preconditioning, solve

  • symbolically represents the preconditioning matrix or process

  • Inside GMRES, we only apply the action of on a vector

  • Right preconditioned matrix free version

PETSc Options

PETSc parameters can either be set on the command line or by using the "petsc_options", "petsc_options_iname", and "petsc_options_value" parameters. Several PETSc parameters that users could frequently use are listed below:

Common stand-alone petsc options

petsc_optionsDescription
-snes_ksp_ewVariable linear solve tolerance – useful for transient solves
-helpShow PETSc options during the solve

Common petsc options and values

petsc_options_inamepetsc_options_valueDescription
-pc_typeiluDefault for serial
bjacobiDefault for parallel with -sub_pc_type ilu
asmAdditive Schwartz with -sub_pc_type ilu
luFull LU, serial only!
gamgPETSc Geometric AMG Preconditioner
hypreHypre, usually used with boomeramg
-sub_pc_typeilu, lu, hypreCan be used with bjacobi or asm
-pc_hypre_typeboomeramgAlgebraic multigrid
-pc_hypre_boomeramg (cont.)"Information Threshold" for AMG process
_strong_threshold0.0 - 1.0(0.7 is auto set for 3D)
-ksp_gmres_restart#Number of Krylov vectors to store

Input Parameters

  • splittingTop-level splitting defining a hierarchical decomposition into subsystems to help the solver.

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Top-level splitting defining a hierarchical decomposition into subsystems to help the solver.

  • time0System time

    Default:0

    C++ Type:double

    Controllable:No

    Description:System time

  • verboseFalseSet to true to print additional information

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Set to true to print additional information

Optional Parameters

  • accept_on_max_fixed_point_iterationFalseTrue to treat reaching the maximum number of fixed point iterations as converged.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:True to treat reaching the maximum number of fixed point iterations as converged.

  • auto_advanceFalseWhether to automatically advance sub-applications regardless of whether their solve converges, for transient executioners only.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to automatically advance sub-applications regardless of whether their solve converges, for transient executioners only.

  • custom_abs_tol1e-50The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on postprocessor defined by the custom_pp residual.

    Default:1e-50

    C++ Type:double

    Controllable:No

    Description:The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on postprocessor defined by the custom_pp residual.

  • custom_ppPostprocessor for custom fixed point convergence check.

    C++ Type:PostprocessorName

    Controllable:No

    Description:Postprocessor for custom fixed point convergence check.

  • custom_rel_tol1e-08The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the postprocessor defined by custom_pp residual.

    Default:1e-08

    C++ Type:double

    Controllable:No

    Description:The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the postprocessor defined by custom_pp residual.

  • direct_pp_valueFalseTrue to use direct postprocessor value (scaled by value on first iteration). False (default) to use difference in postprocessor value between fixed point iterations.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:True to use direct postprocessor value (scaled by value on first iteration). False (default) to use difference in postprocessor value between fixed point iterations.

  • disable_fixed_point_residual_norm_checkFalseDisable the residual norm evaluation thus the three parameters fixed_point_rel_tol, fixed_point_abs_tol and fixed_point_force_norms.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Disable the residual norm evaluation thus the three parameters fixed_point_rel_tol, fixed_point_abs_tol and fixed_point_force_norms.

  • fixed_point_abs_tol1e-50The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

    Default:1e-50

    C++ Type:double

    Controllable:No

    Description:The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

  • fixed_point_algorithmpicardThe fixed point algorithm to converge the sequence of problems.

    Default:picard

    C++ Type:MooseEnum

    Options:picard, secant, steffensen

    Controllable:No

    Description:The fixed point algorithm to converge the sequence of problems.

  • fixed_point_force_normsFalseForce the evaluation of both the TIMESTEP_BEGIN and TIMESTEP_END norms regardless of the existence of active MultiApps with those execute_on flags, default: false.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Force the evaluation of both the TIMESTEP_BEGIN and TIMESTEP_END norms regardless of the existence of active MultiApps with those execute_on flags, default: false.

  • fixed_point_max_its1Specifies the maximum number of fixed point iterations.

    Default:1

    C++ Type:unsigned int

    Controllable:No

    Description:Specifies the maximum number of fixed point iterations.

  • fixed_point_min_its1Specifies the minimum number of fixed point iterations.

    Default:1

    C++ Type:unsigned int

    Controllable:No

    Description:Specifies the minimum number of fixed point iterations.

  • fixed_point_rel_tol1e-08The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

    Default:1e-08

    C++ Type:double

    Controllable:No

    Description:The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

  • relaxation_factor1Fraction of newly computed value to keep.Set between 0 and 2.

    Default:1

    C++ Type:double

    Controllable:No

    Description:Fraction of newly computed value to keep.Set between 0 and 2.

  • transformed_postprocessorsList of main app postprocessors to transform during fixed point iterations

    C++ Type:std::vector<PostprocessorName>

    Controllable:No

    Description:List of main app postprocessors to transform during fixed point iterations

  • transformed_variablesList of main app variables to transform during fixed point iterations

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:List of main app variables to transform during fixed point iterations

Fixed Point Iterations Parameters

  • automatic_scalingFalseWhether to use automatic scaling for the variables.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to use automatic scaling for the variables.

  • compute_scaling_onceTrueWhether the scaling factors should only be computed once at the beginning of the simulation through an extra Jacobian evaluation. If this is set to false, then the scaling factors will be computed during an extra Jacobian evaluation at the beginning of every time step.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Whether the scaling factors should only be computed once at the beginning of the simulation through an extra Jacobian evaluation. If this is set to false, then the scaling factors will be computed during an extra Jacobian evaluation at the beginning of every time step.

  • ignore_variables_for_autoscalingList of variables that do not participate in autoscaling.

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:List of variables that do not participate in autoscaling.

  • off_diagonals_in_auto_scalingFalseWhether to consider off-diagonals when determining automatic scaling factors.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to consider off-diagonals when determining automatic scaling factors.

  • resid_vs_jac_scaling_param0A parameter that indicates the weighting of the residual vs the Jacobian in determining variable scaling parameters. A value of 1 indicates pure residual-based scaling. A value of 0 indicates pure Jacobian-based scaling

    Default:0

    C++ Type:double

    Controllable:No

    Description:A parameter that indicates the weighting of the residual vs the Jacobian in determining variable scaling parameters. A value of 1 indicates pure residual-based scaling. A value of 0 indicates pure Jacobian-based scaling

  • scaling_group_variablesName of variables that are grouped together for determining scale factors. (Multiple groups can be provided, separated by semicolon)

    C++ Type:std::vector<std::vector<std::string>>

    Controllable:No

    Description:Name of variables that are grouped together for determining scale factors. (Multiple groups can be provided, separated by semicolon)

Solver Variable Scaling Parameters

  • contact_line_search_allowed_lambda_cuts2The number of times lambda is allowed to be cut in half in the contact line search. We recommend this number be roughly bounded by 0 <= allowed_lambda_cuts <= 3

    Default:2

    C++ Type:unsigned int

    Controllable:No

    Description:The number of times lambda is allowed to be cut in half in the contact line search. We recommend this number be roughly bounded by 0 <= allowed_lambda_cuts <= 3

  • contact_line_search_ltolThe linear relative tolerance to be used while the contact state is changing between non-linear iterations. We recommend that this tolerance be looser than the standard linear tolerance

    C++ Type:double

    Controllable:No

    Description:The linear relative tolerance to be used while the contact state is changing between non-linear iterations. We recommend that this tolerance be looser than the standard linear tolerance

  • line_searchdefaultSpecifies the line search type (Note: none = basic)

    Default:default

    C++ Type:MooseEnum

    Options:basic, bt, contact, cp, default, l2, none, project, shell

    Controllable:No

    Description:Specifies the line search type (Note: none = basic)

  • line_search_packagepetscThe solver package to use to conduct the line-search

    Default:petsc

    C++ Type:MooseEnum

    Options:petsc, moose

    Controllable:No

    Description:The solver package to use to conduct the line-search

Solver Line Search Parameters

  • control_tagsAdds user-defined labels for accessing object parameters via control logic.

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Adds user-defined labels for accessing object parameters via control logic.

  • enableTrueSet the enabled status of the MooseObject.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Set the enabled status of the MooseObject.

  • outputsVector of output names where you would like to restrict the output of variables(s) associated with this object

    C++ Type:std::vector<OutputName>

    Controllable:No

    Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object

  • skip_exception_checkFalseSpecifies whether or not to skip exception check

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Specifies whether or not to skip exception check

Advanced Parameters

  • l_abs_tol1e-50Linear Absolute Tolerance

    Default:1e-50

    C++ Type:double

    Controllable:No

    Description:Linear Absolute Tolerance

  • l_max_its10000Max Linear Iterations

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:Max Linear Iterations

  • l_tol1e-05Linear Relative Tolerance

    Default:1e-05

    C++ Type:double

    Controllable:No

    Description:Linear Relative Tolerance

  • reuse_preconditionerFalseIf true reuse the previously calculated preconditioner for the linearized system across multiple solves spanning nonlinear iterations and time steps. The preconditioner resets as controlled by reuse_preconditioner_max_linear_its

    Default:False

    C++ Type:bool

    Controllable:No

    Description:If true reuse the previously calculated preconditioner for the linearized system across multiple solves spanning nonlinear iterations and time steps. The preconditioner resets as controlled by reuse_preconditioner_max_linear_its

  • reuse_preconditioner_max_linear_its25Reuse the previously calculated preconditioner for the linear system until the number of linear iterations exceeds this number

    Default:25

    C++ Type:unsigned int

    Controllable:No

    Description:Reuse the previously calculated preconditioner for the linear system until the number of linear iterations exceeds this number

Linear Solver Parameters

  • max_xfem_update4294967295Maximum number of times to update XFEM crack topology in a step due to evolving cracks

    Default:4294967295

    C++ Type:unsigned int

    Controllable:No

    Description:Maximum number of times to update XFEM crack topology in a step due to evolving cracks

  • update_xfem_at_timestep_beginFalseShould XFEM update the mesh at the beginning of the timestep

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Should XFEM update the mesh at the beginning of the timestep

Xfem Fixed Point Iterations Parameters

  • mffd_typewpSpecifies the finite differencing type for Jacobian-free solve types. Note that the default is wp (for Walker and Pernice).

    Default:wp

    C++ Type:MooseEnum

    Options:wp, ds

    Controllable:No

    Description:Specifies the finite differencing type for Jacobian-free solve types. Note that the default is wp (for Walker and Pernice).

  • petsc_optionsSingleton PETSc options

    C++ Type:MultiMooseEnum

    Options:-dm_moose_print_embedding, -dm_view, -ksp_converged_reason, -ksp_gmres_modifiedgramschmidt, -ksp_monitor, -ksp_monitor_snes_lg-snes_ksp_ew, -ksp_snes_ew, -snes_converged_reason, -snes_ksp, -snes_ksp_ew, -snes_linesearch_monitor, -snes_mf, -snes_mf_operator, -snes_monitor, -snes_test_display, -snes_view

    Controllable:No

    Description:Singleton PETSc options

  • petsc_options_inameNames of PETSc name/value pairs

    C++ Type:MultiMooseEnum

    Options:-ksp_atol, -ksp_gmres_restart, -ksp_max_it, -ksp_pc_side, -ksp_rtol, -ksp_type, -mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -snes_atol, -snes_linesearch_type, -snes_ls, -snes_max_it, -snes_rtol, -snes_divergence_tolerance, -snes_type, -sub_ksp_type, -sub_pc_type

    Controllable:No

    Description:Names of PETSc name/value pairs

  • petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname"

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname"

Petsc Parameters

  • n_max_nonlinear_pingpong100The maximum number of times the nonlinear residual can ping pong before requesting halting the current evaluation and requesting timestep cut

    Default:100

    C++ Type:unsigned int

    Controllable:No

    Description:The maximum number of times the nonlinear residual can ping pong before requesting halting the current evaluation and requesting timestep cut

  • nl_abs_div_tol1e+50Nonlinear Absolute Divergence Tolerance. A negative value disables this check.

    Default:1e+50

    C++ Type:double

    Controllable:No

    Description:Nonlinear Absolute Divergence Tolerance. A negative value disables this check.

  • nl_abs_step_tol0Nonlinear Absolute step Tolerance

    Default:0

    C++ Type:double

    Controllable:No

    Description:Nonlinear Absolute step Tolerance

  • nl_abs_tol1e-50Nonlinear Absolute Tolerance

    Default:1e-50

    C++ Type:double

    Controllable:No

    Description:Nonlinear Absolute Tolerance

  • nl_div_tol1e+10Nonlinear Relative Divergence Tolerance. A negative value disables this check.

    Default:1e+10

    C++ Type:double

    Controllable:No

    Description:Nonlinear Relative Divergence Tolerance. A negative value disables this check.

  • nl_forced_its0The Number of Forced Nonlinear Iterations

    Default:0

    C++ Type:unsigned int

    Controllable:No

    Description:The Number of Forced Nonlinear Iterations

  • nl_max_funcs10000Max Nonlinear solver function evaluations

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:Max Nonlinear solver function evaluations

  • nl_max_its50Max Nonlinear Iterations

    Default:50

    C++ Type:unsigned int

    Controllable:No

    Description:Max Nonlinear Iterations

  • nl_rel_step_tol0Nonlinear Relative step Tolerance

    Default:0

    C++ Type:double

    Controllable:No

    Description:Nonlinear Relative step Tolerance

  • nl_rel_tol1e-08Nonlinear Relative Tolerance

    Default:1e-08

    C++ Type:double

    Controllable:No

    Description:Nonlinear Relative Tolerance

  • num_grids1The number of grids to use for a grid sequencing algorithm. This includes the final grid, so num_grids = 1 indicates just one solve in a time-step

    Default:1

    C++ Type:unsigned int

    Controllable:No

    Description:The number of grids to use for a grid sequencing algorithm. This includes the final grid, so num_grids = 1 indicates just one solve in a time-step

  • residual_and_jacobian_togetherFalseWhether to compute the residual and Jacobian together.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to compute the residual and Jacobian together.

  • snesmf_reuse_baseTrueSpecifies whether or not to reuse the base vector for matrix-free calculation

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Specifies whether or not to reuse the base vector for matrix-free calculation

  • solve_typePJFNK: Preconditioned Jacobian-Free Newton Krylov JFNK: Jacobian-Free Newton Krylov NEWTON: Full Newton Solve FD: Use finite differences to compute Jacobian LINEAR: Solving a linear problem

    C++ Type:MooseEnum

    Options:PJFNK, JFNK, NEWTON, FD, LINEAR

    Controllable:No

    Description:PJFNK: Preconditioned Jacobian-Free Newton Krylov JFNK: Jacobian-Free Newton Krylov NEWTON: Full Newton Solve FD: Use finite differences to compute Jacobian LINEAR: Solving a linear problem

Nonlinear Solver Parameters

    Restart Parameters

    Input Files

    Child Objects