SteadyAndAdjoint

Executioner for evaluating steady-state simulations and their adjoint.

Overview

This executioner can be used to solve a steady-state forward problem with its adjoint. The forward solve is performed the same way as in Steady, but with the "solve_type" set to NEWTON. This performs the nonlinear iteration in the form:

where is the number of iterations it took to converge the problem. The adjoint problem is then defined as:

where is the adjoint solution, is the residual of the adjoint system with , and is the transpose of the forward system's Jacobian evaluated with the converged forward solution. The adjoint system is basically a linearized and homogenized version of the forward problem with it's own definition of sources.

In order to accurately define the adjoint system, the fully consistent Jacobian must be evaluated. As such, the computeQpJacobian routines in the forward problem kernels must be accurately defined or Automatic Differentiation must be used. Consider using the Jacobian debugger to ensure the Jacobian is computed accurately.

Example Input File Syntax

Solving Adjoint Problems

The first step is to add an adjoint nonlinear system using the "nl_sys_names" parameter in the Problem input block. It is convenient to define the forward system as nl0.

[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
  nl_sys_names = 'nl0 adjoint'
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/self_adjoint.i)

Next we need to add an adjoint variable for each forward variable, which is associated with the adjoint system:

Single adjoint variable

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [u]
  []
  [u_adjoint]
    solver_sys = adjoint
  []
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/self_adjoint.i)

Multiple adjoint variables

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [u]
  []
  [v]
  []
  [u_adjoint]
    solver_sys = adjoint
  []
  [v_adjoint]
    solver_sys = adjoint
  []
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/multi_variable.i)

Array adjoint variables

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [u]
    components = 2
  []
  [u_adjoint]
    components = 2
    solver_sys = adjoint
  []
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/array_variable.i)

Next we add kernels and BCs associated with the forward and adjoint variables. Only source-like kernels should be added to the adjoint variables like BodyForce, ConstantPointSource, or NeumannBC.

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [diff]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
  []
  [src]
    type = BodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../kernels/BodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 1
  []
  [src_adjoint]
    type = BodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../kernels/BodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u_adjoint
    value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 10
  []
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [dirichlet]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top right'
    value<<<{"description": "Value of the BC"}>>> = 1
  []
  [neumann]
    type = NeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h$, where $h$ is a constant, controllable value.", "href": "../bcs/NeumannBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left bottom'
    value<<<{"description": "For a Laplacian problem, the value of the gradient dotted with the normals on the boundary."}>>> = 1
  []
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/nonhomogeneous_bc.i)

For nonlinear, problems one should use AD Kernels, BCs, and Materials.

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [diff]
    type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../kernels/ADMatDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    diffusivity<<<{"description": "The diffusivity value or material property"}>>> = D
  []
  [src]
    type = ADBodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../kernels/BodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 1
  []

  [src_adjoint]
    type = ADBodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../kernels/BodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u_adjoint
    value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 1
  []
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [dirichlet]
    type = ADDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/ADDirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top right'
    value<<<{"description": "Value of the BC"}>>> = 0
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [diffc]
    type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../materials/ParsedMaterial.html"}>>>
    property_name<<<{"description": "Name of the parsed material property"}>>> = D
    expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '0.1 + 5 * u'
    coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'u'
  []
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/nonlinear_diffusion.i)

Finally, we will add this executioner and set the forward/adjoint system tags. Note that the tolerance for the adjoint system solve is set solely by linear solver parameters like "l_tol", "l_abs_tol", and "l_max_its".

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = SteadyAndAdjoint
  forward_system = nl0
  adjoint_system = adjoint

  nl_rel_tol = 1e-12
  l_tol = 1e-12
[]
(modules/optimization/test/tests/executioners/steady_and_adjoint/nonhomogeneous_bc.i)

Utilization in Gradient-Based Optimization

Utilizing this executioner for gradient-based optimization is quite powerful since the adjoint used to compute the gradient is automatically assembled with this executioner, given that the forward problem Jacobian can be fully constructed.

To include the source for the adjoint problem, the ReporterPointSource can be used to add the simulation misfit from the forward solve, which is calculated in OptimizationData.

[Reporters<<<{"href": "../../syntax/Reporters/index.html"}>>>]
  [measurement_locations]
    type = OptimizationData<<<{"description": "Reporter to hold measurement and simulation data for optimization problems", "href": "../reporters/OptimizationData.html"}>>>
    objective_name<<<{"description": "Name of reporter value defining the objective."}>>> = objective_value
    variable<<<{"description": "Vector of variable names to sample at measurement points."}>>> = forwardT
  []
  [params]
    type = ConstantReporter<<<{"description": "Reporter with constant values to be accessed by other objects, can be modified using transfers.", "href": "../reporters/ConstantReporter.html"}>>>
    real_vector_names<<<{"description": "Names for each vector of reals value."}>>> = 'heat_source'
    real_vector_values<<<{"description": "Values for vectors of reals."}>>> = '0' # Dummy
  []
[]

[DiracKernels<<<{"href": "../../syntax/DiracKernels/index.html"}>>>]
  [pt]
    type = ReporterPointSource<<<{"description": "Apply a point load defined by Reporter.", "href": "../dirackernels/ReporterPointSource.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = adjointT
    x_coord_name<<<{"description": "reporter x-coordinate name.  This uses the reporter syntax <reporter>/<name>."}>>> = measurement_locations/measurement_xcoord
    y_coord_name<<<{"description": "reporter y-coordinate name.  This uses the reporter syntax <reporter>/<name>."}>>> = measurement_locations/measurement_ycoord
    z_coord_name<<<{"description": "reporter z-coordinate name.  This uses the reporter syntax <reporter>/<name>."}>>> = measurement_locations/measurement_zcoord
    value_name<<<{"description": "reporter value name.  This uses the reporter syntax <reporter>/<name>."}>>> = measurement_locations/misfit_values
  []
[]
(modules/optimization/test/tests/optimizationreporter/nonlinear_material/forward_and_adjoint.i)

The gradient can then be computed using an inner-product vector-postprocessor like ElementOptimizationSourceFunctionInnerProduct. Note that these vector-postprocessors must be executed on ADJOINT_TIMESTEP_END which occurs after the adjoint system is solved.

[VectorPostprocessors<<<{"href": "../../syntax/VectorPostprocessors/index.html"}>>>]
  [gradient_vpp]
    type = ElementOptimizationSourceFunctionInnerProduct<<<{"description": "Computes the inner product of variable with parameterized source function for optimization gradient computation.", "href": "../vectorpostprocessors/ElementOptimizationSourceFunctionInnerProduct.html"}>>>
    function<<<{"description": "Optimization function."}>>> = volumetric_heat_func
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = adjointT
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = ADJOINT_TIMESTEP_END
  []
[]
(modules/optimization/test/tests/optimizationreporter/nonlinear_material/forward_and_adjoint.i)

The driving optimize app will thus have only have a single FullSolveMultiApp which then transfers both the simulation data (for the objective evaluation) and the inner products (for the gradient evaluation).

[MultiApps<<<{"href": "../../syntax/MultiApps/index.html"}>>>]
  [forward]
    type = FullSolveMultiApp<<<{"description": "Performs a complete simulation during each execution.", "href": "../multiapps/FullSolveMultiApp.html"}>>>
    input_files<<<{"description": "The input file for each App.  If this parameter only contains one input file it will be used for all of the Apps.  When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = forward_and_adjoint.i
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = FORWARD
  []
[]

[Transfers<<<{"href": "../../syntax/Transfers/index.html"}>>>]
  [to_forward]
    type = MultiAppReporterTransfer<<<{"description": "Transfers reporter data between two applications.", "href": "../transfers/MultiAppReporterTransfer.html"}>>>
    to_multi_app<<<{"description": "The name of the MultiApp to transfer the data to"}>>> = forward
    from_reporters<<<{"description": "List of the reporter names (object_name/value_name) to transfer the value from."}>>> = 'main/measurement_xcoord
                      main/measurement_ycoord
                      main/measurement_zcoord
                      main/measurement_time
                      main/measurement_values
                      OptimizationReporter/heat_source'
    to_reporters<<<{"description": "List of the reporter names (object_name/value_name) to transfer the value to."}>>> = 'measurement_locations/measurement_xcoord
                    measurement_locations/measurement_ycoord
                    measurement_locations/measurement_zcoord
                    measurement_locations/measurement_time
                    measurement_locations/measurement_values
                    params/heat_source'
  []
  [from_forward]
    type = MultiAppReporterTransfer<<<{"description": "Transfers reporter data between two applications.", "href": "../transfers/MultiAppReporterTransfer.html"}>>>
    from_multi_app<<<{"description": "The name of the MultiApp to receive data from"}>>> = forward
    from_reporters<<<{"description": "List of the reporter names (object_name/value_name) to transfer the value from."}>>> = 'measurement_locations/objective_value
                      gradient_vpp/inner_product'
    to_reporters<<<{"description": "List of the reporter names (object_name/value_name) to transfer the value to."}>>> = 'OptimizationReporter/objective_value
                    OptimizationReporter/grad_heat_source'
  []
[]
(modules/optimization/test/tests/optimizationreporter/nonlinear_material/main.i)

Input Parameters

  • adjoint_systemName of the system representing the adjoint problem.

    C++ Type:NonlinearSystemName

    Controllable:No

    Description:Name of the system representing the adjoint problem.

  • forward_systemName of the nonlinear system representing the forward problem. Multiple and linear solver systems are not currently supported.

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

    Controllable:No

    Description:Name of the nonlinear system representing the forward problem. Multiple and linear solver systems are not currently supported.

Required Parameters

  • time0System time

    Default:0

    C++ Type:double

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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.

  • multiapp_fixed_point_convergenceName of the Convergence object to use to assess convergence of the MultiApp fixed point solve. If not provided, a default Convergence will be constructed internally from the executioner parameters.

    C++ Type:ConvergenceName

    Controllable:No

    Description:Name of the Convergence object to use to assess convergence of the MultiApp fixed point solve. If not provided, a default Convergence will be constructed internally from the executioner parameters.

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

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    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>

    Unit:(no unit assumed)

    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_once1 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. Vector entries correspond to each nonlinear system.

    Default:1

    C++ Type:std::vector<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. Vector entries correspond to each nonlinear system.

  • ignore_variables_for_autoscalingList of variables that do not participate in autoscaling. Vector entries correspond to each nonlinear system.

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

    Controllable:No

    Description:List of variables that do not participate in autoscaling. Vector entries correspond to each nonlinear system.

  • off_diagonals_in_auto_scaling0 Whether to consider off-diagonals when determining automatic scaling factors. Vector entries correspond to each nonlinear system.

    Default:0

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

    Controllable:No

    Description:Whether to consider off-diagonals when determining automatic scaling factors. Vector entries correspond to each nonlinear system.

  • resid_vs_jac_scaling_param0 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. Vector entries correspond to each nonlinear system.

    Default:0

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

    Unit:(no unit assumed)

    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. Vector entries correspond to each nonlinear system.

  • scaling_group_variablesName of variables that are grouped together for determining scale factors. (Multiple groups can be provided, separated by semicolon). Vector entries correspond to each nonlinear system.

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

    Controllable:No

    Description:Name of variables that are grouped together for determining scale factors. (Multiple groups can be provided, separated by semicolon). Vector entries correspond to each nonlinear system.

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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

  • linear_convergenceName of the Convergence object(s) to use to assess convergence of the linear system(s) solve. If not provided, the linear solver tolerance parameters are used

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

    Controllable:No

    Description:Name of the Convergence object(s) to use to assess convergence of the linear system(s) solve. If not provided, the linear solver tolerance parameters are used

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

    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 for transient simulations

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

    Default:1e+50

    C++ Type:double

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    Controllable:No

    Description:Nonlinear Absolute step Tolerance

  • nl_abs_tol1e-50Nonlinear Absolute Tolerance

    Default:1e-50

    C++ Type:double

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    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

    Unit:(no unit assumed)

    Controllable:No

    Description:Nonlinear Relative step Tolerance

  • nl_rel_tol1e-08Nonlinear Relative Tolerance

    Default:1e-08

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:Nonlinear Relative Tolerance

  • nonlinear_convergenceName of the Convergence object(s) to use to assess convergence of the nonlinear system(s) solve. If not provided, the default Convergence associated with the Problem will be constructed internally.

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

    Controllable:No

    Description:Name of the Convergence object(s) to use to assess convergence of the nonlinear system(s) solve. If not provided, the default Convergence associated with the Problem will be constructed internally.

  • 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_together0 Whether to compute the residual and Jacobian together. Vector entries correspond to each nonlinear system.

    Default:0

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

    Controllable:No

    Description:Whether to compute the residual and Jacobian together. Vector entries correspond to each nonlinear system.

  • 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

  • splittingTop-level splitting defining a hierarchical decomposition into subsystems to help the solver. Outer-vector of this vector-of-vector parameter correspond to each nonlinear system.

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

    Controllable:No

    Description:Top-level splitting defining a hierarchical decomposition into subsystems to help the solver. Outer-vector of this vector-of-vector parameter correspond to each nonlinear system.

  • use_pre_SMO_residualFalseCompute the pre-SMO residual norm and use it in the relative convergence check. The pre-SMO residual is computed at the begining of the time step before solution-modifying objects are executed. Solution-modifying objects include preset BCs, constraints, predictors, etc.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Compute the pre-SMO residual norm and use it in the relative convergence check. The pre-SMO residual is computed at the begining of the time step before solution-modifying objects are executed. Solution-modifying objects include preset BCs, constraints, predictors, etc.

Nonlinear 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_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

    Controllable:No

    Description:Singleton PETSc options

  • petsc_options_inameNames of PETSc name/value pairs

    C++ Type:MultiMooseEnum

    Options:-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, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_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

  • multi_system_fixed_pointFalseWhether to perform fixed point (Picard) iterations between the nonlinear systems.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to perform fixed point (Picard) iterations between the nonlinear systems.

  • multi_system_fixed_point_convergenceConvergence object to determine the convergence of the multi-system fixed point iteration. If unspecified, defaults to checking that every system is converged (based on their own convergence criterion)

    C++ Type:ConvergenceName

    Controllable:No

    Description:Convergence object to determine the convergence of the multi-system fixed point iteration. If unspecified, defaults to checking that every system is converged (based on their own convergence criterion)

Multiple Solver System Parameters

    Restart Parameters

    Input Files