[petsc-users] replacing Schur complement in pcfieldsplit

Klaij, Christiaan C.Klaij at marin.nl
Thu Jun 21 02:17:49 CDT 2012


> >
> > In order to implement SIMPLE-type preconditioners for the
> > incompressible Navier-Stokes equations (Elman e.a. JCP 227, 2008)
> > using the PCFieldSplit framework, it looks like I need to replace
> > inv(A00) by some cheap approximation
> >
> > 1) in the Schur complement
> >
>
> When you have a Schur FS, the '0' solver is this approximation.
>
>
> > 2) in the L and/or U of the LDU factorization
> >
>
> You can use whatever PC you want for the solver mentioned above.
>
>
> > 3) while keeping A00 in the D
> >
>
> I think what you want here is -pc_fieldsplit_real_diagonal.

Let me get this straight. Looking at the full LDU factorization
of the block matrix. Citing from the manual:

  For the Schur complement preconditioner if
  J = ( A00 A01 )
      ( A10 A11 )

  the preconditioner using full factorization is
      ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
      ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )

Clearly inv(A00) occurs four times, right? In L and in U and
twice in D. Now if I somehow overrule the '0' solver with my
approximation and use -pc_fieldsplit_real_diagonal, the effect
would be that inv(A00) is replaced in L, in U and in S but not in
the 00-block of D?

And what's the function name corresponding to
-pc_fieldsplit_real_diagonal?

>
> I would be really nice to make a PETSc example that did SIMPLE. I was
> going to do this, but if you do it first, I will put it in.

Working on it.


dr. ir. Christiaan Klaij
CFD Researcher
Research & Development
E mailto:C.Klaij at marin.nl
T +31 317 49 33 44

MARIN
2, Haagsteeg, P.O. Box 28, 6700 AA Wageningen, The Netherlands
T +31 317 49 39 11, F +31 317 49 32 45, I www.marin.nl



More information about the petsc-users mailing list