[petsc-users] Convergence problem in solving a symmetric positive definite matrix in a CFD code

Matthew Knepley knepley at gmail.com
Wed Dec 4 09:19:56 CST 2013


On Wed, Dec 4, 2013 at 8:59 AM, Xiao, Jianjun (IKET)
<jianjun.xiao at kit.edu>wrote:

> Hello,
>
> I am using Petsc to solve the linear equation in a CFD code. The matrix is
> symmetric positive definite.
>
> Please find my input and output below.
>
> 1. When I used KSPCR solver, ||r(i)||/||b|| is inf. What is the reason?
> Does it mean ||b|| is zero? When I used the KSPLGMRES solver, it seems
> ||r(i)||/||b|| is OK. However, it seems the calculated results are not
> right.
>

It looks that way. You can check it with VecNorm() before calling
KSPSolve().

   Matt


> 2. I am not sure if I set the solver and matrix properly. Did I miss
> something?
>
> Thank you.
>
>
> ****************************************
> INPUT:
>
>        CALL KSPCreate(PETSC_COMM_WORLD,solver,ierr)
>        CALL DMSetMatType(da_ksp,MATMPISBAIJ,ierr)
>        CALL DMCreateMatrix(da_ksp,mat,ierr)
>        CALL
> MatMPISBAIJSetPreallocation(mat,1,4,PETSC_NULL_INTEGER,1,PETSC_NULL_INTEGER,ierr)
>        CALL MatZeroEntries(mat,ierr)
>        CALL MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY,ierr)
>        CALL MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY,ierr)
>
>        CALL KSPSetOperators(solver,mat,mat,SAME_NONZERO_PATTERN,ierr)
>
> !      CALL KSPSetType(solver,KSPCR,ierr)
>        CALL KSPSetType(solver,KSPLGMRES,ierr)
>
>        epsi = 1.0-e5
>        CALL KSPSetTolerances(solver,epsi,PETSC_DEFAULT_DOUBLE_PRECISION,&
>            & PETSC_DEFAULT_DOUBLE_PRECISION,itmax,ierr)
>
>        CALL KSPGetPC(solver,gfprec,ierr)
>        CALL PCSetType(prec,PCBJACOBI,ierr)
>
>        CALL KSPMonitorSet(solver,KSPMonitorTrueResidualNorm,&
>         PETSC_NULL_OBJECT,PETSC_NULL_FUNCTION,ierr )
>
>        CALL KSPSetInitialGuessNonzero(solver,PETSC_TRUE,ierr)
>        CALL KSPSolve(solver,bvec,xsol,ierr)
>
>        CALL KSPGetIterationNumber(solver,iter,ierr)
>        CALL KSPGetResidualNorm(solver,dmax,ierr)
>        CALL KSPView(solver,PETSC_VIEWER_STDOUT_WORLD,ierr)
>
> **************************************
> OUTPUT CALL KSPSetType(solver,KSPCR,ierr):
>
>   0 KSP preconditioned resid norm 2.226482634319e+05 true resid norm
> 1.204978940624e+07 ||r(i)||/||b||            inf
>   1 KSP preconditioned resid norm 1.684243557244e+05 true resid norm
> 6.742321430949e+06 ||r(i)||/||b||            inf
>   2 KSP preconditioned resid norm 1.039386033131e+05 true resid norm
> 5.094347016880e+06 ||r(i)||/||b||            inf
>   3 KSP preconditioned resid norm 3.767761162917e+04 true resid norm
> 2.085014289432e+06 ||r(i)||/||b||            inf
>   4 KSP preconditioned resid norm 2.220316358489e+04 true resid norm
> 1.039841616110e+06 ||r(i)||/||b||            inf
>   5 KSP preconditioned resid norm 1.009108756815e+04 true resid norm
> 6.764592620620e+05 ||r(i)||/||b||            inf
>   6 KSP preconditioned resid norm 7.266143334386e+03 true resid norm
> 4.713756053613e+05 ||r(i)||/||b||            inf
>   7 KSP preconditioned resid norm 4.925270026573e+03 true resid norm
> 3.276759177651e+05 ||r(i)||/||b||            inf
>   8 KSP preconditioned resid norm 2.595039666791e+03 true resid norm
> 1.774916597474e+05 ||r(i)||/||b||            inf
>   9 KSP preconditioned resid norm 1.970388137453e+03 true resid norm
> 1.449811653036e+05 ||r(i)||/||b||            inf
>  10 KSP preconditioned resid norm 1.455914234388e+03 true resid norm
> 7.916294162841e+04 ||r(i)||/||b||            inf
>  11 KSP preconditioned resid norm 8.335194818556e+02 true resid norm
> 4.530953608250e+04 ||r(i)||/||b||            inf
>  12 KSP preconditioned resid norm 3.022320555777e+02 true resid norm
> 1.728551635678e+04 ||r(i)||/||b||            inf
>  13 KSP preconditioned resid norm 7.190336024797e+01 true resid norm
> 4.186842086105e+03 ||r(i)||/||b||            inf
>  14 KSP preconditioned resid norm 1.718291655675e+01 true resid norm
> 1.089751055004e+03 ||r(i)||/||b||            inf
>  15 KSP preconditioned resid norm 1.150683059424e+01 true resid norm
> 8.672405273471e+02 ||r(i)||/||b||            inf
>  16 KSP preconditioned resid norm 8.663479440949e+00 true resid norm
> 5.776737380768e+02 ||r(i)||/||b||            inf
>  17 KSP preconditioned resid norm 5.282161990683e+00 true resid norm
> 2.977735906695e+02 ||r(i)||/||b||            inf
>  18 KSP preconditioned resid norm 3.802629315725e+00 true resid norm
> 2.789114564993e+02 ||r(i)||/||b||            inf
>  19 KSP preconditioned resid norm 1.722575171383e+00 true resid norm
> 1.051323829526e+02 ||r(i)||/||b||            inf
> KSP Object: 1 MPI processes
>   type: cr
>   maximum iterations=1000
>   tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>   left preconditioning
>   using nonzero initial guess
>   using PRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>   type: bjacobi
>     block Jacobi: number of blocks = 1
>     Local solve is same for all blocks, in the following KSP and PC
> objects:
>     KSP Object:    (sub_)     1 MPI processes
>       type: preonly
>       maximum iterations=10000, initial guess is zero
>       tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>       left preconditioning
>       using NONE norm type for convergence test
>     PC Object:    (sub_)     1 MPI processes
>       type: icc
>         0 levels of fill
>         tolerance for zero pivot 2.22045e-14
>         using Manteuffel shift [POSITIVE_DEFINITE]
>         matrix ordering: natural
>         factor fill ratio given 1, needed 1
>           Factored matrix follows:
>             Mat Object:             1 MPI processes
>               type: seqsbaij
>               rows=27, cols=27
>               package used to perform factorization: petsc
>               total: nonzeros=81, allocated nonzeros=81
>               total number of mallocs used during MatSetValues calls =0
>                   block size is 1
>       linear system matrix = precond matrix:
>       Mat Object:       1 MPI processes
>         type: seqsbaij
>         rows=27, cols=27
>         total: nonzeros=81, allocated nonzeros=108
>         total number of mallocs used during MatSetValues calls =0
>             block size is 1
>   linear system matrix = precond matrix:
>   Mat Object:   1 MPI processes
>     type: mpisbaij
>     rows=27, cols=27
>     total: nonzeros=81, allocated nonzeros=135
>     total number of mallocs used during MatSetValues calls =0
>         block size is 1
>
> OUTPUT: CALL KSPSetType(solver,KSPLGMRES,ierr)
>
>   0 KSP preconditioned resid norm 2.362537325084e+04 true resid norm
> 1.138584383312e+06 ||r(i)||/||b|| 1.000000000000e+00
>   1 KSP preconditioned resid norm 8.501213683423e+03 true resid norm
> 3.655528853686e+05 ||r(i)||/||b|| 3.210591070162e-01
>   2 KSP preconditioned resid norm 5.487567253725e+03 true resid norm
> 3.005741194777e+05 ||r(i)||/||b|| 2.639893220769e-01
>   3 KSP preconditioned resid norm 2.470452880657e+03 true resid norm
> 1.545469272201e+05 ||r(i)||/||b|| 1.357360328187e-01
>   4 KSP preconditioned resid norm 1.750803325456e+03 true resid norm
> 1.182312309352e+05 ||r(i)||/||b|| 1.038405520646e-01
>   5 KSP preconditioned resid norm 1.123492053552e+03 true resid norm
> 6.754319630701e+04 ||r(i)||/||b|| 5.932208213726e-02
>   6 KSP preconditioned resid norm 5.150241959277e+02 true resid norm
> 3.689413898730e+04 ||r(i)||/||b|| 3.240351749775e-02
>   7 KSP preconditioned resid norm 4.182894544871e+02 true resid norm
> 3.052196222024e+04 ||r(i)||/||b|| 2.680693909701e-02
>   8 KSP preconditioned resid norm 2.520093155629e+02 true resid norm
> 1.880976788356e+04 ||r(i)||/||b|| 1.652031079932e-02
>   9 KSP preconditioned resid norm 1.186491314806e+02 true resid norm
> 6.797080217853e+03 ||r(i)||/||b|| 5.969764136483e-03
>  10 KSP preconditioned resid norm 5.851092821372e+01 true resid norm
> 2.973659280245e+03 ||r(i)||/||b|| 2.611716201127e-03
>  11 KSP preconditioned resid norm 1.669909189055e+01 true resid norm
> 5.658829814125e+02 ||r(i)||/||b|| 4.970057465277e-04
>  12 KSP preconditioned resid norm 3.090594692756e+00 true resid norm
> 2.161527454147e+02 ||r(i)||/||b|| 1.898434130864e-04
>  13 KSP preconditioned resid norm 2.164618839184e+00 true resid norm
> 1.620745991834e+02 ||r(i)||/||b|| 1.423474637093e-04
>  14 KSP preconditioned resid norm 1.291593952428e+00 true resid norm
> 9.095542547366e+01 ||r(i)||/||b|| 7.988465923722e-05
>  15 KSP preconditioned resid norm 6.100583411632e-01 true resid norm
> 4.021646656091e+01 ||r(i)||/||b|| 3.532146334551e-05
>  16 KSP preconditioned resid norm 2.723496807925e-01 true resid norm
> 1.676660466866e+01 ||r(i)||/||b|| 1.472583403954e-05
>  17 KSP preconditioned resid norm 1.377718471538e-01 true resid norm
> 8.551854245272e+00 ||r(i)||/||b|| 7.510953400221e-06
> KSP Object: 1 MPI processes
>   type: lgmres
>     GMRES: restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>     GMRES: happy breakdown tolerance 1e-30
>     LGMRES: aug. dimension=2
>     LGMRES: number of matvecs=170
>   maximum iterations=1000
>   tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>   left preconditioning
>   using nonzero initial guess
>   using PRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>   type: bjacobi
>     block Jacobi: number of blocks = 1
>     Local solve is same for all blocks, in the following KSP and PC
> objects:
>     KSP Object:    (sub_)     1 MPI processes
>       type: preonly
>       maximum iterations=10000, initial guess is zero
>       tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>       left preconditioning
>       using NONE norm type for convergence test
>     PC Object:    (sub_)     1 MPI processes
>       type: icc
>         0 levels of fill
>         tolerance for zero pivot 2.22045e-14
>         using Manteuffel shift [POSITIVE_DEFINITE]
>         matrix ordering: natural
>         factor fill ratio given 1, needed 1
>           Factored matrix follows:
>             Mat Object:             1 MPI processes
>               type: seqsbaij
>               rows=27, cols=27
>               package used to perform factorization: petsc
>               total: nonzeros=81, allocated nonzeros=81
>               total number of mallocs used during MatSetValues calls =0
>                   block size is 1
>       linear system matrix = precond matrix:
>       Mat Object:       1 MPI processes
>         type: seqsbaij
>         rows=27, cols=27
>         total: nonzeros=81, allocated nonzeros=108
>         total number of mallocs used during MatSetValues calls =0
>             block size is 1
>   linear system matrix = precond matrix:
>   Mat Object:   1 MPI processes
>     type: mpisbaij
>     rows=27, cols=27
>     total: nonzeros=81, allocated nonzeros=135
>     total number of mallocs used during MatSetValues calls =0
>



-- 
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/20131204/4271255f/attachment.html>


More information about the petsc-users mailing list