<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Dec 13, 2013 at 4:16 AM, 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">Hello Jed and Matt,<br>
<br>
Thanks for your reply.<br>
<br>
I design a small 3D case with 10 cells in each direction. Please find the binary files generated with PETSc-dev in the attachment. The matrix should be symmetric and positive definite.<br>
<br>
The calculation did not converge at the first time step with PCR+JACOBI. With PCG+JACOBI, it seems it converges in this case. I have very limited experience in using PETSc. In the past, I felt more flexible to change the combination of linear solvers and preconditioners. I could get divergence anyhow, and the main difference was the speed. But now, it seems only PCG+JACOBI works fine for my case. Many other combinations, such as PCR+JACOBI, PCR+BJACOBI,  PCG+BJACOBI, divergent at the first time step or take much longer time. I think there must be something wrong.<br>
</blockquote><div><br></div><div>Can you see that the above paragraph gives NO IDEA what your problem is?</div><div><br></div><div>1) You completely failed to do the one thing I asked you to do. That is not the way to get help.</div>
<div><br></div><div>     *** Run the same solver with 3.3 and dev. Does it give the same answer ***</div><div><br></div><div>    Until you do this, why would I believe anything has changed.</div><div><br></div><div>2) Iterative solvers fail to converge all the time. That is routine, and the reason for preconditioning. Preconditioners</div>
<div>     are always tied to the specific problem. Thus telling us that something did not converge without outlining the</div><div>     problem, or providing the output I asked for (-ksp_monitor_true_residual -ksp_converged_reason) tells us</div>
<div>     nothing.</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">
Below is how I set up the matrix and the solver.<br>
<br>
       CALL DMSetMatType(da_ksp,MATMPISBAIJ,ierr)<br>
       CALL DMCreateMatrix(da_ksp,gmat,ierr)<br>
       CALL DMSetMatrixPreallocateOnly(da_ksp,PETSC_TRUE,ierr)<br>
       CALL MatMPISBAIJSetPreallocation(gmat,1,4,PETSC_NULL_INTEGER,3,PETSC_NULL_INTEGER,ierr)<br>
       CALL MatZeroEntries(gmat,ierr)<br>
<br>
       CALL KSPCreate(PETSC_COMM_WORLD,gsolver,ierr)<br>
!       CALL KSPSetType(gsolver,KSPCR,ierr)<br>
       CALL KSPSetType(gsolver,KSPCG,ierr)<br>
<br>
       CALL KSPGetPC(gsolver,gprec,ierr)<br>
       CALL PCSetType(gprec,PCBJACOBI,ierr)<br>
!       CALL PCSetType(gprec,PCBJACOBI,ierr)<br>
<br>
********************************** matsetvalues<br>
       ...<br>
           IF (k < kbar-1 .AND. j < jbar-1 .AND. i < ibar-1) THEN ! d-x-y-z<br>
<br>
             valpm(1)     = value_ad_gfpp(m)<br>
             colpm(MatStencil_i,1) = i<br>
             colpm(MatStencil_j,1) = j<br>
             colpm(MatStencil_k,1) = k<br>
<br>
             valpm(2)     = value_ax_gfpp(m)<br>
             colpm(MatStencil_i,2) = i+1<br>
             colpm(MatStencil_j,2) = j<br>
             colpm(MatStencil_k,2) = k<br>
<br>
             valpm(3)     = value_ay_gfpp(m)<br>
             colpm(MatStencil_i,3) = i<br>
             colpm(MatStencil_j,3) = j+1<br>
             colpm(MatStencil_k,3) = k<br>
<br>
             valpm(4)     = value_az_gfpp(m)<br>
             colpm(MatStencil_i,4) = i<br>
             colpm(MatStencil_j,4) = j<br>
             colpm(MatStencil_k,4) = k+1<br>
<br>
             CALL MatSetValuesBlockedStencil(gmat,1,rowpm,4,colpm,valpm,INSERT_VALUES,ierr)<br>
       ...<br>
<br>
       ...<br>
**********************************<br>
       CALL MatAssemblyBegin(gmat,MAT_FINAL_ASSEMBLY,ierr)<br>
       CALL MatAssemblyEnd(gmat,MAT_FINAL_ASSEMBLY,ierr)<br>
<br>
       CALL KSPSetOperators(gsolver,gmat,gmat,SAME_NONZERO_PATTERN,ierr)<br>
<br>
       epsi_rtol = 1.0e-50_dble<br>
       epsi_atol = 1.0e-5_dble<br>
<br>
       CALL KSPSetTolerances(gsolver,epsi_rtol, epsi_atol,&<br>
           & PETSC_DEFAULT_DOUBLE_PRECISION,itmax,ierr)<br>
<br>
       CALL KSPSetNormType(gsolver,KSP_NORM_PRECONDITIONED,ierr)<br>
<br>
       CALL KSPSolve(gsolver,bvec0_gfpp,xsol0_gfpp,ierr)<br>
<br>
Thank you for your help.<br>
<br>
JJ<br>
<br>
<br>
________________________________________<br>
From: Jed Brown [<a href="mailto:jedbrown@mcs.anl.gov">jedbrown@mcs.anl.gov</a>]<br>
Sent: Friday, December 13, 2013 8:34 AM<br>
To: Xiao, Jianjun (IKET); Matthew Knepley<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>
<div class="HOEnZb"><div class="h5"><br>
"Xiao, Jianjun (IKET)" <<a href="mailto:jianjun.xiao@kit.edu">jianjun.xiao@kit.edu</a>> writes:<br>
<br>
> Hi Matt,<br>
><br>
> I am using PETSc-dev because only this version supports<br>
> MatSetValuesBlockedStencil in Fortran. I tried many combinations of<br>
> KSP and PC. None of them worked in my case. The same thing happened to<br>
> PETSc-3.4.<br>
><br>
> Then I changed back to 3.3-P7 and 3.2-P7 without using<br>
> MatSetValuesBlockedStencil. I tried KSPCR+PCBJACOBI and many others,<br>
> they worked fine.<br>
<br>
We're going to need exact options and output for the cases you tried.<br>
You should be able to use the same algorithm with both versions (with<br>
suitable options).  And if you have a modest-sized problem, you can run<br>
with -ksp_view_mat and send binaryoutput.  We'll use<br>
src/ksp/ksp/examples/tutorials/ex10.c to reproduce the convergence you<br>
see.<br>
</div></div></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>