petsc-3.9.4 2018-09-11
Report Typos and Errors

KSPTSIRM

Implements the two-stage iteration with least-squares residual minimization method.

Options Database Keys

-ksp_ksp_type <solver> - the type of the inner solver (GMRES or any of its variants for instance)
-ksp_pc_type <preconditioner> - the type of the preconditioner applied to the inner solver
-ksp_ksp_max_it <maxits> - the maximum number of inner iterations (iterations of the inner solver)
-ksp_ksp_rtol <tol> - sets the relative convergence tolerance of the inner solver
-ksp_tsirm_cgls <number> - if 1 use CGLS solver in the minimization step, otherwise use LSQR solver
-ksp_tsirm_max_it_ls <maxits> - the maximum number of iterations for the least-squares minimization solver
-ksp_tsirm_tol_ls <tol> - sets the convergence tolerance of the least-squares minimization solver
-ksp_tsirm_size_ls <size> - the number of residuals for the least-squares minimization step

Notes: TSIRM is a new two-stage iteration method for solving large sparse linear systems of the form Ax=b. The main idea behind this new method is the use a least-squares residual minimization to improve the convergence of Krylov based iterative methods, typically those of GMRES variants. The principle of TSIRM algorithm is to build an outer iteration over a Krylov method, called inner solver, and to frequently store the current residual computed by the given Krylov method in a matrix of residuals S. After a few outer iterations, a least-squares minimization step is applied on the matrix composed by the saved residuals, in order to compute a better solution and to make new iterations if required. The GMRES method , or any of its variants, can potentially be used as inner solver. The minimization step consists in solving the least-squares problem min||b-ASa|| to find 'a' which minimizes the residuals (b-AS). The minimization step is performed using two solvers of linear least-squares problems: CGLS or LSQR. A new solution x with a minimal residual is computed with x=Sa.

References

1 R. Couturier, L. Ziane Khodja, and C. Guyeux. TSIRM: A Two -Stage Iteration with least-squares Residual Minimization algorithm to solve large sparse linear systems. In PDSEC 2015, 16th IEEE Int. Workshop on Parallel and Distributed Scientific and Engineering Computing (in conjunction with IPDPS 2015), Hyderabad, India, 2015.

Contributed by: Lilia Ziane Khodja

See Also

KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()

Level

advanced

Location

src/ksp/ksp/impls/tsirm/tsirm.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages