- topsplitEntrance to splits, the top split will specify how splits will go.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Entrance to splits, the top split will specify how splits will go.
FSP
Preconditioner designed to map onto PETSc's PCFieldSplit.
Overview
The FieldSplitPreconditioner
allows for custom preconditioning for each nonlinear variable in the numerical system. One or more variables may be targeted in a subsolve that will only consider part of the numerical system. The preconditioning defined for these subsolves is used for the relevant block(s) in the global numerical system.
A FSP
may for example be used for block-diagonal preconditioning by setting full=false
and no off-diagonal variable couplings. Numerical systems considering only a single variable are then preconditioned individually. This is the default preconditioner for the PJFNK
solves. See the Executioner documentation for more information on the default preconditioner.
More information about field split preconditioning may be found in the PETSc manual.
Example input syntax
In this example, the preconditioning is performed by solving individual problems for each variables, as described in the comments in the snippet. The solution for these subsolves is used to perform the Schur decomposition preconditioning of the main numerical system.
[Preconditioning]
active = 'FSP'
[./FSP]
type = FSP
# It is the starting point of splitting
topsplit = 'uv' # 'uv' should match the following block name
[./uv]
splitting = 'u v' # 'u' and 'v' are the names of subsolvers
# Generally speaking, there are four types of splitting we could choose
# <additive,multiplicative,symmetric_multiplicative,schur>
splitting_type = additive
# An approximate solution to the original system
# | A_uu A_uv | | u | _ |f_u|
# | 0 A_vv | | v | - |f_v|
# is obtained by solving the following subsystems
# A_uu u = f_u and A_vv v = f_v
# If splitting type is specified as schur, we may also want to set more options to
# control how schur works using PETSc options
# petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition'
# petsc_options_value = 'full selfp'
[../]
[./u]
vars = 'u'
# PETSc options for this subsolver
# A prefix will be applied, so just put the options for this subsolver only
petsc_options_iname = '-pc_type -ksp_type'
petsc_options_value = ' hypre preonly'
[../]
[./v]
vars = 'v'
# PETSc options for this subsolver
petsc_options_iname = '-pc_type -ksp_type'
petsc_options_value = ' hypre preonly'
[../]
[../]
[]
(moose/test/tests/preconditioners/fsp/fsp_test.i)An example of setting the "off_diag_row" and "off_diag_column" parameters to create a custom coupling matrix may be found in the PBP documentation.
Input Parameters
- fullTrueSet to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination.
- ksp_normunpreconditionedSets the norm that is used for convergence testing
Default:unpreconditioned
C++ Type:MooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Sets the norm that is used for convergence testing
- nl_sysThe nonlinear system whose linearization this preconditioner should be applied to.
C++ Type:NonlinearSystemName
Unit:(no unit assumed)
Controllable:No
Description:The nonlinear system whose linearization this preconditioner should be applied to.
- off_diag_columnThe variable names for the off-diagonal columns you want to add into the matrix; they will be associated with an off-diagonal row from the same position in off_diag_row.
C++ Type:std::vector<NonlinearVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The variable names for the off-diagonal columns you want to add into the matrix; they will be associated with an off-diagonal row from the same position in off_diag_row.
- off_diag_rowThe variable names for the off-diagonal rows you want to add into the matrix; they will be associated with an off-diagonal column from the same position in off_diag_column.
C++ Type:std::vector<NonlinearVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The variable names for the off-diagonal rows you want to add into the matrix; they will be associated with an off-diagonal column from the same position in off_diag_column.
- pc_sidedefaultPreconditioning side
Default:default
C++ Type:MooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Preconditioning side
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
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
Unit:(no unit assumed)
Controllable:No
Description:Set the enabled status of the MooseObject.
Advanced 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
Unit:(no unit assumed)
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
Unit:(no unit assumed)
Controllable:No
Description:Singleton PETSc options
- petsc_options_inameNames of PETSc name/value pairs
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
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>
Unit:(no unit assumed)
Controllable:No
Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname"
- 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
Unit:(no unit assumed)
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