[petsc-users] MINRES and PCFIELDSPLIT

Matthew Knepley knepley at gmail.com
Mon Jan 20 09:01:45 CST 2014


On Mon, Jan 20, 2014 at 8:38 AM, Adriano Côrtes <adrimacortes at gmail.com>wrote:

> Hi all,
>
> I'm still trying to make MINRES and PCFIELDSPLIT working together in
> my code, but I am not understanding how things are going.
> Ok, first of all, I'm trying to solve a Stokes problem with an inf-sup
> stable element with all boundary of Dirichlet type, that is, my system
> is
>
>        [ A   B^t ]
> M = [            ]
>        [ B    0   ]
>
> and its kernel is the column vector [zeros(ndof_vel,1) ;
> ones(ndof_press,1) ] (using Matlab notation)
>
> Lets say some facts that I'm sure of, since I damped the matrices
> using MatView to a Matlab file and inspected them:
>
> a. M is symmetric for sure and singular, with the kernel mentioned above;
> b. The block A is SPD (all eigenvalues are positive);
> c. S = - B * A^{-1} * B^t negative semi-definite (with only one zero
> eigenvalue and all others negative); that is -S is positive
> semi-definite;
>
> Now I will talk about my tests using MINRES and PCFIELDSPLIT to solve
> this problem:
>
> 1. MINRES without preconditioner converges (very slowly) to any
> specified tolerance;
>
> 2. MINRES without preconditioner AND using KSPSetNullspace to set the
> null space of the matrix as [zeros(ndof_vel,1) ; ones(ndof_press,1) ]
> converges to any specified tolerance, but the number of iterations are
> almost the same of the test above;
>
> 3. MINRES with PCFIELDSPLIT, using SPLIT = SCHUR, and SCHURFACT =
> DIAG, and letting Petsc managing the other options with default values
> (GMRES for the KPS and ILU(0) for PC). In this case I think Petsc is
> using
>
>        [  A    0 ]
> P  = [           ]
>        [  0   -S ]
>
> that in my case is symmetric positive semi-definite (by just one zero
> eigenvalue), but approximating with ksp(A) and ksp(-S). Before someone
> as I am also setting  KSPSetNullspace as above for the MINRES
> iterations, plus passing the option
> -fieldsplit_p_ksp_constant_null_space for the Schur block, and what I
> obtain (with -ksp_rtol 1e-7) is
>
>     0 KSP preconditioned resid norm 4.632440509011e-01 true resid norm
> 2.987970930100e-02 ||r(i)||/||b|| 1.000000000000e+00
>     1 KSP preconditioned resid norm 1.484245122855e-01 true resid norm
> 3.127586084436e-03 ||r(i)||/||b|| 1.046725740512e-01
>     2 KSP preconditioned resid norm 1.474500219584e-01 true resid norm
> 3.303438991154e-03 ||r(i)||/||b|| 1.105579360855e-01
>     3 KSP preconditioned resid norm 3.685525320417e-06 true resid norm
> 2.519862718263e-07 ||r(i)||/||b|| 8.433357543339e-06
>     4 KSP preconditioned resid norm 3.665488598091e-06 true resid norm
> 2.739588290212e-07 ||r(i)||/||b|| 9.168724710853e-06
>     5 KSP preconditioned resid norm 9.845910914885e-07 true resid norm
> 2.763504619849e-07 ||r(i)||/||b|| 9.248766753419e-06
>     6 KSP preconditioned resid norm 9.519073659560e-08 true resid norm
> 1.490709049851e-07 ||r(i)||/||b|| 4.989034648343e-06
>     7 KSP preconditioned resid norm 9.490881783501e-08 true resid norm
> 1.489003256545e-07 ||r(i)||/||b|| 4.983325779864e-06
>   Linear solve did not converge due to DIVERGED_INDEFINITE_MAT iterations 8
>
> So I ask: why this last line? What is computed inside Petsc that
> forces it to be printed? My guess is that o preconditioned norm used.
> Am I right?
>

The problem here is the definition of S. In order for the theory to work, S
needs to be
accurate, so you should always start by solving A exactly,

  -fieldsplit_u_ksp_type preonly
  -fieldsplit_y_pc_type lu

If after you see 3 iterate (I think) convergence, then you can back off
this solver,
but you likely need to solve A inside S accurately, which you can do using

  -fieldsplit_p_inner_ksp_type preonly
  -fieldsplit_p_inner_pc_type lu

See http://people.cs.uchicago.edu/~knepley/presentations/SIAMCSE13.pdf. You
can then backoff LU and see the change in convergence. Eventually, if the
solve
is inaccurate enough, you can destroy the definiteness.

   Matt

Then I decide to change some options for PCFIELDSPLIT, passing the options
>
> #Velocity
> -fieldsplit_u_ksp_type cg
> -fieldsplit_u_pc_type bjacobi
> -fieldsplit_u_ksp_rtol 1e-03
>
> #Pressure
> -fieldsplit_p_ksp_type cg
> -fieldsplit_p_pc_type bjacobi
> -fieldsplit_p_ksp_rtol 1e-03
> -fieldsplit_p_ksp_constant_null_space
>
> Then I obtained the -ksp_monitor (again with -ksp_rtol 1e-7):
>
>     0 KSP preconditioned resid norm 4.632522584339e-01 true resid norm
> 2.987970930100e-02 ||r(i)||/||b|| 1.000000000000e+00
>     1 KSP preconditioned resid norm 1.485092747354e-01 true resid norm
> 3.133707103392e-03 ||r(i)||/||b|| 1.048774294229e-01
>     2 KSP preconditioned resid norm 1.475544621957e-01 true resid norm
> 3.307719603712e-03 ||r(i)||/||b|| 1.107011976051e-01
>     3 KSP preconditioned resid norm 3.647747431146e-04 true resid norm
> 3.795610417928e-05 ||r(i)||/||b|| 1.270296969657e-03
>     4 KSP preconditioned resid norm 3.614914908198e-04 true resid norm
> 4.114945598888e-05 ||r(i)||/||b|| 1.377170559939e-03
>     5 KSP preconditioned resid norm 1.289144645113e-04 true resid norm
> 4.067177049739e-05 ||r(i)||/||b|| 1.361183607500e-03
>     6 KSP preconditioned resid norm 5.143993331403e-05 true resid norm
> 2.794092581519e-05 ||r(i)||/||b|| 9.351137099001e-04
>     7 KSP preconditioned resid norm 2.361302329725e-05 true resid norm
> 2.715477815463e-05 ||r(i)||/||b|| 9.088032912597e-04
>     8 KSP preconditioned resid norm 2.352126784416e-05 true resid norm
> 2.732009930616e-05 ||r(i)||/||b|| 9.143361814852e-04
>   Linear solve did not converge due to DIVERGED_INDEFINITE_MAT iterations 9
>
> And again the same line is printed. Can somebody explain it for me?
> Thanks.
>
> Adriano Cortes
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140120/47b7aed6/attachment-0001.html>


More information about the petsc-users mailing list