[petsc-users] replacing Schur complement in pcfieldsplit

Klaij, Christiaan C.Klaij at marin.nl
Fri Jun 22 02:46:06 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  )
> >
>
> Yes.
>
>
> > 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
> >
>
> Yes
>
>
> > 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?
> >
>
> No. What this says is that we should use the action of the
> actual matrix rather than the preconditioner matrix in the solver.

Odd. Don't you *always* need the action of the actual matrix (and
of the preconditioner) in a Krylov subspace method?

>
> I now think I have a better idea what you want, but it would be
> helpful if you wrote it out completely in linear algebra notation, as
> above. Right now, we use the same KSP for all 4 places above.
> Using different KSPs would require a small code change, which I
> can make if you give me a better idea what you want.

Maybe there is a mistake in the manual, shouldn't it be

 ( I   -A01 ksp(A00) )
 ( 0         I       )

in the factorization above instead of -A10 ksp(A00)? SIMPLE-type
preconditioners are usually written as:

 ( I   -A01 dA00^(-1) ) ( A00 0 )^(-1)
 ( 0         I        ) ( A10 S )

with S = A11 - A10 dA00^(-1) A01, where dA00^(-1) is the inverse
of the diagonal of A00. Therefore it only requires one solve for
A00 and one solve for S.

>
>
> > And what's the function name corresponding to
> > -pc_fieldsplit_real_diagonal?
> >
>
> We have not put one in yet.

Please let me know when you do.

>
>    Thanks,
>
>       Matt


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