<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Thu, Dec 12, 2013 at 3:09 PM, Xiao, Jianjun (IKET) <span dir="ltr"><<a href="mailto:jianjun.xiao@kit.edu" target="_blank">jianjun.xiao@kit.edu</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi Matt,<br>
<br>
I am using PETSc-dev because only this version supports MatSetValuesBlockedStencil in Fortran. I tried  many combinations of KSP and PC. None of them worked in my case. The same thing happened to PETSc-3.4.<br>
<br>
Then I changed back to 3.3-P7 and 3.2-P7 without using MatSetValuesBlockedStencil. I tried KSPCR+PCBJACOBI and many others,  they worked fine.<br>
<br>
I read the changes of each version, and still could not find the reason. Did I miss something? Thank you.<br></blockquote><div><br></div><div>1) "Did not work" is an insufficient description of the problem. If it did not converge, give the output of -ksp_monitor_true_residual -ksp_converged_reason</div>
<div>      for each run</div><div><br></div><div>2) You changed two things at once. That makes it impossible to debug. Use petsc-dev without MatSetValuesBlockedStencil(). Do you</div><div>     reproduce the results from 3.3? If so, I believe your blocks may be ordered incorrectly.</div>
<div><br></div><div>   Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
________________________________________<br>
From: Matthew Knepley [<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>]<br>
Sent: Wednesday, December 04, 2013 4:19 PM<br>
To: Xiao, Jianjun (IKET)<br>
Cc: <a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a><br>
Subject: Re: [petsc-users] Convergence problem in solving a symmetric positive definite matrix in a CFD code<br>
<br>
On Wed, Dec 4, 2013 at 8:59 AM, Xiao, Jianjun (IKET) <<a href="mailto:jianjun.xiao@kit.edu">jianjun.xiao@kit.edu</a><mailto:<a href="mailto:jianjun.xiao@kit.edu">jianjun.xiao@kit.edu</a>>> wrote:<br>
Hello,<br>
<br>
I am using Petsc to solve the linear equation in a CFD code. The matrix is symmetric positive definite.<br>
<br>
Please find my input and output below.<br>
<br>
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.<br>

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