KSPFETIDP#
The FETI-DP method [FLL+01]
Options Database Keys#
-ksp_fetidp_fullyredundant
- use a fully redundant set of Lagrange multipliers-ksp_fetidp_saddlepoint
- activates support for saddle point problems, see [TL15]-ksp_fetidp_saddlepoint_flip
- see note below-ksp_fetidp_pressure_field <- 1> - activates support for saddle point problems, and identifies the pressure field id. If this information is not provided, the pressure field is detected by using
MatFindZeroDiagonals()
.-ksp_fetidp_pressure_all
- if false, uses the interface pressures, as described in [2]. If true, uses the entire pressure field.
Notes#
The matrix for the KSP
must be of type MATIS
.
Usually, an incompressible Stokes problem is written as
| A B^T | | v | = | f |
| B 0 | | p | = | g |
with B representing \( -\int_\Omega \nabla \cdot u q \). If -ksp_fetidp_saddlepoint_flip is true, the code assumes that the user provides it as
| A B^T | | v | = | f |
|-B 0 | | p | = |-g |
The FETI-DP linear system (automatically generated constructing an internal PCBDDC
object) is solved using an internal KSP
object.
Options for the inner KSP
and for the customization of the PCBDDC
object can be specified at command line by using the prefixes -fetidp_
and -fetidp_bddc_
. E.g.,
-fetidp_ksp_type gmres -fetidp_bddc_pc_bddc_symmetric false
will use KSPGMRES
for the solution of the linear system on the Lagrange multipliers, generated using a non-symmetric PCBDDC
.
For saddle point problems with continuous pressures, the preconditioned operator for the pressure solver can be specified with KSPFETIDPSetPressureOperator()
.
Alternatively, the pressure operator is extracted from the precondioned matrix (if it is different from the linear solver matrix).
If none of the above, an identity matrix will be created; the user then needs to scale it through a Richardson solver.
Options for the pressure solver can be prefixed with -fetidp_fielsplit_p_
, E.g.
-fetidp_fielsplit_p_ksp_type preonly -fetidp_fielsplit_p_pc_type lu -fetidp_fielsplit_p_pc_factor_mat_solver_type mumps
In order to use the deluxe version of FETI-DP, you must customize the inner PCBDDC
operator with -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat and use
non-redundant multipliers, i.e. -ksp_fetidp_fullyredundant false
. Options for the scaling solver are prefixed by -fetidp_bddelta_
, E.g.
-fetidp_bddelta_pc_factor_mat_solver_type mumps -fetidp_bddelta_pc_type lu
Some of the basic options such as the maximum number of iterations and tolerances are automatically passed from this KSP
to the inner KSP
that actually performs the iterations.
The converged reason and number of iterations computed are passed from the inner KSP
to this KSP
at the end of the solution.
Developer Note#
Even though this method does not directly use any norms, the user is allowed to set the KSPNormType
to any value.
This is so users do not have to change KSPNormType
options when they switch from other KSP
methods to this one.
References#
- FLL+01
Charbel Farhat, Michel Lesoinne, Patrick LeTallec, Kendall Pierson, and Daniel Rixen. FETI-DP: a dual–primal unified FETI method—part I: a faster alternative to the two-level FETI method. International journal for numerical methods in engineering, 50(7):1523–1544, 2001.
- TL15
Xuemin Tu and Jing Li. A FETI-DP type domain decomposition algorithm for three-dimensional incompressible Stokes equations. SIAM Journal on Numerical Analysis, 53(2):720–742, 2015.
See Also#
KSP: Linear System Solvers, MATIS
, PCBDDC
, KSPFETIDPSetInnerBDDC()
, KSPFETIDPGetInnerBDDC()
, KSPFETIDPGetInnerKSP()
Level#
Advanced
Location#
src/ksp/ksp/impls/fetidp/fetidp.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages