[petsc-users] Solve Linear System with Field Split Preconditioner
晓峰 何
tlanyan at hotmail.com
Mon Sep 26 04:20:03 CDT 2022
Are there other approaches to solve this kind of systems in PETSc except for field-split methods?
Thanks,
Xiaofeng
On Sep 26, 2022, at 14:13, Jed Brown <jed at jedbrown.org<mailto:jed at jedbrown.org>> wrote:
This is the joy of factorization field-split methods. The actual Schur complement is dense, so we represent it implicitly. A common strategy is to assemble the mass matrix and drop it in the 11 block of the Pmat. You can check out some examples in the repository for incompressible flow (Stokes problems). The LSC (least squares commutator) is another option. You'll likely find that lumping diag(A00)^{-1} works poorly because the resulting operator behaves like a Laplacian rather than like a mass matrix.
晓峰 何 <tlanyan at hotmail.com<mailto:tlanyan at hotmail.com>> writes:
If assigned a preconditioner to A11 with this cmd options:
-fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type ilu -fieldsplit_1_ksp_type gmres -fieldsplit_1_pc_type ilu
Then I got this error:
"Could not locate a solver type for factorization type ILU and matrix type schurcomplement"
How could I specify a preconditioner for A11?
BR,
Xiaofeng
On Sep 26, 2022, at 11:02, 晓峰 何 <tlanyan at hotmail.com<mailto:tlanyan at hotmail.com><mailto:tlanyan at hotmail.com>> wrote:
-fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type ilu -fieldsplit_1_ksp_type gmres -fieldsplit_1_pc_type none
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20220926/851b5e8c/attachment-0001.html>
More information about the petsc-users
mailing list