[petsc-users] usefulness of fieldsplit_type schur and fieldsplit_schur_type user

Jed Brown jed at jedbrown.org
Wed Feb 15 22:58:00 CST 2017


David Nolte <dnolte at dim.uchile.cl> writes:

> Dear all,
>
> I am trying to use PETSc to solve a steady Stokes problem, discretized
> with stable FEM (P2P1) in 2D and 3D. I have been playing around with
> different preconditioners to get the hang of things. For my 2D case, for
> example, the following setup works well, where fieldsplit_0 is
> identified with the velocity dofs and fieldsplit_1 corresponds to the
> pressure dofs, and the user defined Schur approximation is the pressure
> mass matrix:
>
>            -ksp_view
>            -ksp_converged_reason
>            -ksp_type fgmres
>            -ksp_rtol 1.0e-8
>            -ksp_monitor
>
>             -pc_type fieldsplit
>             -pc_fieldsplit_type schur
>             -pc_fieldsplit_schur_fact_type upper
>             -pc_fieldsplit_schur_precondition user
>
>             -fieldsplit_0_ksp_type richardson
>             -fieldsplit_0_ksp_max_it 1
>             -fieldsplit_0_pc_type hypre
>             -fieldsplit_0_pc_hypre_type boomeramg
>
>             -fieldsplit_1_ksp_type preonly
>             -fieldsplit_1_pc_type jacobi
>
> However, now I wonder what is the difference between
>     -pc_fieldsplit_type schur
>     -pc_fieldsplit_schur_fact_type diag
>     -pc_fieldsplit_schur_precondition user
> and
>     -pc_fieldsplit_type additive
> where the preconditioner matrix would be P = [A, 0; 0, Ŝ] (Ŝ is the same
> Schur approximation as above)? And likewise, between "schur_fact_type
> lower" and "fieldsplit_type multiplicative", with P = [A, 0; B, Ŝ].
> (Except for the generality/flexibility of the Schur approach.)
>
> Is there any advantage of using the Schur approach over the additive
> oder multiplicative fieldsplits if I specify the approx. Schur
> complement Ŝ? Any reason to prefer one method over the other? The Schur
> factorization causes a (very) slightly longer computation time.

Depends if you might want to iterate on the Schur complement.  With
Schur, the operator for the inner solver becomes the Schur complement S
= D - C A^{-1} B where A^{-1} is an iterative solve.  This is Schur
complement reduction.  If you use preonly, this representation isn't
used.

With additive or multiplicative, that operator is just the block of the
outer operator (D) which is nothing like the Schur complement.  So the
schur option is better for debugging convergence or if you want to
iterate on the Schur complement.
-------------- next part --------------
A non-text attachment was scrubbed...
Name: signature.asc
Type: application/pgp-signature
Size: 832 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20170215/b49363f8/attachment.pgp>


More information about the petsc-users mailing list