<html><head><meta http-equiv="Content-Type" content="text/html; charset=utf-8"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div><div class=""><div class="">Before you do any investigation I would run with one SNES solve </div><div class=""><br class=""></div><div class="">-snes_max_it 1 -snes_view </div><div class=""><br class=""></div><div class="">It will print out exactly what solver configuration you are using.</div></div><div class=""><br class=""></div><div class="">------</div><div class=""><br class=""></div><div class="">[0]PETSC ERROR: KSPSolve has not converged, reason DIVERGED_DTOL</div><div class="">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" class="">https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.</div><div class="">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020 </div><div class="">[0]PETSC ERROR: ./stokeTutorial on a arch-linux2-c-debug named a2aa8f1c96aa by Unknown Fri Oct  9 13:43:28 2020</div><div class="">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda</div><div class="">[0]PETSC ERROR: #1 KSPSolve() line 832 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c</div><div class="">[0]PETSC ERROR: #2 PCApply_FieldSplit_Schur() line 1189 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c</div><div class=""><br class=""></div><div class="">  Ok, one of the inner solvers went crazy and diverged, that is the norm of the residual for that inner solver exploded. </div><div class=""><br class=""></div><div class="">  Base on the line number </div><div class=""><br class=""></div><div class="">[0]PETSC ERROR: #2 PCApply_FieldSplit_Schur() line 1189 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c </div><div class=""><br class=""></div><div class="">  you can look at that file and see which of the inner solvers failed. From the info from -snes_view you will know what the KSP and PC is for that inner solve and the KSP options prefix. With that you can run the failing case with the addition option </div><div class=""><br class=""></div><div class="">  -xxx_ksp_monitor_true_residual </div><div class=""><br class=""></div><div class="">and watch the inner solve explode. </div><div class=""><br class=""></div><div class="">This inner solver behaving badly can also explain the need for -ksp_type fgmres.  Normally PCFIELDSPLIT is a linear operator and so you can use -ksp_type gmres but there is some issue with the inner solver. Could it possible have a null space, that is be singular? Do you provide your own custom inner solver or just select from the options database? Does -pc_type lu make everything work fine?</div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class="">  </div><div><br class=""><blockquote type="cite" class=""><div class="">On Oct 9, 2020, at 12:53 AM, Karl Yang <<a href="mailto:y.juntao@hotmail.com" class="">y.juntao@hotmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><meta http-equiv="Content-Type" content="text/html; charset=utf-8" class=""><div class="">Hi, Barry, </div><br class=""><div class="">Thanks for your reply. Yes, I should have used fgmres. But after switching to fgmres I'm still facing the same convergence issue. </div><br class=""><div class="">Seems like the reason is due to DIVERGED_PC_FAILED. But I simply used FD jacobian, and fieldsplitPC.  I am a bit lost on whether I made some mistakes somewhere in the FormFunction or I did not setup the solver correctly. </div><br class=""><div class="">///////code///////</div><div class="">  SNESSetFunction(snes, r, FormFunctionStatic, this);</div><div class="">   // SNESSetJacobian(snes, J, J, FormJacobianStatic, this);</div><div class="">     SNESSetJacobian(snes, J, J, SNESComputeJacobianDefault, this);  </div><div class="">      SNESMonitorSet(snes, MySNESMonitor, NULL, NULL);</div><br class=""><div class=""> SNESGetKSP(snes, &ksp);</div><div class="">   KSPGetPC(ksp, &pc);</div><div class="">       PCSetType(pc, PCFIELDSPLIT);</div><div class="">  PCFieldSplitSetDetectSaddlePoint(pc, PETSC_TRUE);</div><div class="">     PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_SELF, NULL);</div><div class="">      KSPMonitorSet(ksp, MyKSPMonitor, NULL, 0);</div><div class="">    KSPSetTolerances(ksp, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);</div><div class="">     SNESSetFromOptions(snes);</div><div class="">//////end/////////</div><br class=""><div class="">Output from SNES/KSP solver</div><div class="">################# step 1 #################</div><div class="">iter = 0, SNES Function norm 0.0430713</div><div class="">iteration 0 KSP Residual norm 4.307133784528e-02 </div><div class="">    0 KSP unpreconditioned resid norm 4.307133784528e-02 true resid norm 4.307133784528e-02 ||r(i)||/||b|| 1.000000000000e+00</div><div class="">iteration 1 KSP Residual norm 4.451434065870e-07 </div><div class="">    1 KSP unpreconditioned resid norm 4.451434065870e-07 true resid norm 4.451434065902e-07 ||r(i)||/||b|| 1.033502623460e-05</div><div class="">iteration 2 KSP Residual norm 1.079756105012e-12 </div><div class="">    2 KSP unpreconditioned resid norm 1.079756105012e-12 true resid norm 1.079754870815e-12 ||r(i)||/||b|| 2.506898844643e-11</div><div class="">  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div class="">iter = 1, SNES Function norm 2.40846e-05</div><div class="">iteration 0 KSP Residual norm 2.408462930023e-05 </div><div class="">    0 KSP unpreconditioned resid norm 2.408462930023e-05 true resid norm 2.408462930023e-05 ||r(i)||/||b|| 1.000000000000e+00</div><div class="">iteration 1 KSP Residual norm 1.096958085415e-11 </div><div class="">    1 KSP unpreconditioned resid norm 1.096958085415e-11 true resid norm 1.096958085425e-11 ||r(i)||/||b|| 4.554598170270e-07</div><div class="">iteration 2 KSP Residual norm 5.909523288165e-16 </div><div class="">    2 KSP unpreconditioned resid norm 5.909523288165e-16 true resid norm 5.909519599233e-16 ||r(i)||/||b|| 2.453647729249e-11</div><div class="">  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div class="">iter = 2, SNES Function norm 1.19684e-14</div><div class="">################# step 2 #################</div><div class="">iter = 0, SNES Function norm 0.00391662</div><div class="">iteration 0 KSP Residual norm 3.916615614134e-03 </div><div class="">    0 KSP unpreconditioned resid norm 3.916615614134e-03 true resid norm 3.916615614134e-03 ||r(i)||/||b|| 1.000000000000e+00</div><div class="">iteration 1 KSP Residual norm 4.068800385009e-08 </div><div class="">    1 KSP unpreconditioned resid norm 4.068800385009e-08 true resid norm 4.068800384986e-08 ||r(i)||/||b|| 1.038856192653e-05</div><div class="">iteration 2 KSP Residual norm 8.427513055511e-14 </div><div class="">    2 KSP unpreconditioned resid norm 8.427513055511e-14 true resid norm 8.427497502034e-14 ||r(i)||/||b|| 2.151729537007e-11</div><div class="">  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div class="">iter = 1, SNES Function norm 1.99152e-07</div><div class="">iteration 0 KSP Residual norm 1.991523558528e-07 </div><div class="">    0 KSP unpreconditioned resid norm 1.991523558528e-07 true resid norm 1.991523558528e-07 ||r(i)||/||b|| 1.000000000000e+00</div><div class="">iteration 1 KSP Residual norm 1.413505562549e-13 </div><div class="">    1 KSP unpreconditioned resid norm 1.413505562549e-13 true resid norm 1.413505562550e-13 ||r(i)||/||b|| 7.097609046588e-07</div><div class="">iteration 2 KSP Residual norm 5.165934822520e-18 </div><div class="">    2 KSP unpreconditioned resid norm 5.165934822520e-18 true resid norm 5.165932973227e-18 ||r(i)||/||b|| 2.593960262787e-11</div><div class="">  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div class="">iter = 2, SNES Function norm 1.69561e-16</div><div class="">################# step 3 #################</div><div class="">iter = 0, SNES Function norm 0.00035615</div><div class="">iteration 0 KSP Residual norm 3.561504844171e-04 </div><div class="">    0 KSP unpreconditioned resid norm 3.561504844171e-04 true resid norm 3.561504844171e-04 ||r(i)||/||b|| 1.000000000000e+00</div><div class="">iteration 1 KSP Residual norm 3.701591890269e-09 </div><div class="">    1 KSP unpreconditioned resid norm 3.701591890269e-09 true resid norm 3.701591890274e-09 ||r(i)||/||b|| 1.039333667153e-05</div><div class="">iteration 2 KSP Residual norm 7.832821034843e-15 </div><div class="">    2 KSP unpreconditioned resid norm 7.832821034843e-15 true resid norm 7.832856926692e-15 ||r(i)||/||b|| 2.199311041093e-11</div><div class="">  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div class="">iter = 1, SNES Function norm 1.64671e-09</div><div class="">iteration 0 KSP Residual norm 1.646709543241e-09 </div><div class="">    0 KSP unpreconditioned resid norm 1.646709543241e-09 true resid norm 1.646709543241e-09 ||r(i)||/||b|| 1.000000000000e+00</div><div class="">iteration 1 KSP Residual norm 1.043230469512e-15 </div><div class="">    1 KSP unpreconditioned resid norm 1.043230469512e-15 true resid norm 1.043230469512e-15 ||r(i)||/||b|| 6.335242749968e-07</div><div class="">iteration 1 KSP Residual norm 0.000000000000e+00 </div><div class="">    1 KSP unpreconditioned resid norm 0.000000000000e+00 true resid norm           -nan ||r(i)||/||b||           -nan</div><div class="">  Linear solve did not converge due to DIVERGED_PC_FAILED iterations 1</div><div class="">                 PC_FAILED due to SUBPC_ERROR </div><br class=""><br class=""><div class="">More information from -ksp_error_if_not_converged -info</div><br class=""><div class="">[0] KSPConvergedDefault(): Linear solver has converged. Residual norm 3.303168180659e-07 is less than relative tolerance 1.000000000000e-05 times initial right hand side norm 7.795816360977e-02 at iteration 12</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] KSPConvergedDefault(): Linear solver has converged. Residual norm 2.227610512466e+00 is less than relative tolerance 1.000000000000e-05 times initial right hand side norm 5.453050347652e+05 at iteration 12</div><div class="">[0] PCSetUp(): Leaving PC with identical preconditioner since operator is unchanged</div><div class="">[0] KSPConvergedDefault(): Linear solver is diverging. Initial right hand size norm 9.501675075823e-01, current residual norm 4.894880836662e+04 at iteration 210</div><div class="">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------</div><div class="">[0]PETSC ERROR:   </div><div class="">[0]PETSC ERROR: KSPSolve has not converged, reason DIVERGED_DTOL</div><div class="">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" class="">https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.</div><div class="">[0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020 </div><div class="">[0]PETSC ERROR: ./stokeTutorial on a arch-linux2-c-debug named a2aa8f1c96aa by Unknown Fri Oct  9 13:43:28 2020</div><div class="">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda</div><div class="">[0]PETSC ERROR: #1 KSPSolve() line 832 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c</div><div class="">[0]PETSC ERROR: #2 PCApply_FieldSplit_Schur() line 1189 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c</div><div class="">[0]PETSC ERROR: #3 PCApply() line 444 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c</div><div class="">[0]PETSC ERROR: #4 KSP_PCApply() line 281 in /usr/local/petsc/petsc-3.12.5/include/petsc/private/kspimpl.h</div><div class="">[0]PETSC ERROR: #5 KSPFGMRESCycle() line 166 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/fgmres/fgmres.c</div><div class="">[0]PETSC ERROR: #6 KSPSolve_FGMRES() line 291 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/fgmres/fgmres.c</div><div class="">[0]PETSC ERROR: #7 KSPSolve() line 760 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c</div><div class="">[0]PETSC ERROR: #8 SNESSolve_NEWTONLS() line 225 in /usr/local/petsc/petsc-3.12.5/src/snes/impls/ls/ls.c</div><div class="">[0]PETSC ERROR: #9 SNESSolve() line 4482 in /usr/local/petsc/petsc-3.12.5/src/snes/interface/snes.c</div><div class="">[0] PetscCommDuplicate(): Using internal PETSc communicator 1140850688 -2080374783</div><div class="">SNES Object: 1 MPI processes</div><div class="">  type: newtonls</div><div class="">  maximum iterations=50, maximum function evaluations=10000</div><div class="">  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08</div><div class="">  total number of linear solver iterations=2</div><div class="">  total number of function evaluations=1322</div><div class="">  norm schedule ALWAYS</div><div class="">  Jacobian is built using finite differences one column at a time</div><div class="">  SNESLineSearch Object: 1 MPI processes</div><div class="">    type: bt</div><div class="">      interpolation: cubic</div><div class="">      alpha=1.000000e-04</div><div class="">    maxstep=1.000000e+08, minlambda=1.000000e-12</div><div class="">    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08</div><div class="">    maximum iterations=40</div><div class="">  KSP Object: 1 MPI processes</div><div class="">    type: fgmres</div><div class="">      restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div class="">      happy breakdown tolerance 1e-30</div><div class="">    maximum iterations=10000, initial guess is zero</div><div class="">    tolerances:  relative=1e-08, absolute=1e-50, divergence=10000.</div><div class="">    right preconditioning</div><div class="">    using UNPRECONDITIONED norm type for convergence test</div><div class="">  PC Object: 1 MPI processes</div><div class="">    type: fieldsplit</div><div class="">      FieldSplit with Schur preconditioner, blocksize = 1, factorization FULL</div><div class="">      Preconditioner for the Schur complement formed from S itself</div><div class="">      Split info:</div><div class="">      Split number 0 Defined by IS</div><div class="">      Split number 1 Defined by IS</div><div class="">      KSP solver for A00 block</div><div class="">        KSP Object: (fieldsplit_0_) 1 MPI processes</div><div class="">          type: gmres</div><div class="">            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div class="">            happy breakdown tolerance 1e-30</div><div class="">          maximum iterations=10000, initial guess is zero</div><div class="">          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.</div><div class="">          left preconditioning</div><div class="">          using PRECONDITIONED norm type for convergence test</div><div class="">        PC Object: (fieldsplit_0_) 1 MPI processes</div><div class="">          type: ilu</div><div class="">            out-of-place factorization</div><div class="">            0 levels of fill</div><div class="">            tolerance for zero pivot 2.22045e-14</div><div class="">            matrix ordering: natural</div><div class="">            factor fill ratio given 1., needed 1.</div><div class="">              Factored matrix follows:</div><div class="">                Mat Object: 1 MPI processes</div><div class="">                  type: seqaij</div><div class="">                  rows=512, cols=512</div><div class="">                  package used to perform factorization: petsc</div><div class="">                  total: nonzeros=9213, allocated nonzeros=9213</div><div class="">                  total number of mallocs used during MatSetValues calls=0</div><div class="">                    not using I-node routines</div><div class="">          linear system matrix = precond matrix:</div><div class="">          Mat Object: 1 MPI processes</div><div class="">            type: seqaij</div><div class="">            rows=512, cols=512</div><div class="">            total: nonzeros=9213, allocated nonzeros=9213</div><div class="">            total number of mallocs used during MatSetValues calls=0</div><div class="">              not using I-node routines</div><div class="">      KSP solver for S = A11 - A10 inv(A00) A01 </div><div class="">        KSP Object: (fieldsplit_1_) 1 MPI processes</div><div class="">          type: gmres</div><div class="">            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div class="">            happy breakdown tolerance 1e-30</div><div class="">          maximum iterations=10000, initial guess is zero</div><div class="">          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.</div><div class="">          left preconditioning</div><div class="">          using PRECONDITIONED norm type for convergence test</div><div class="">        PC Object: (fieldsplit_1_) 1 MPI processes</div><div class="">          type: none</div><div class="">          linear system matrix = precond matrix:</div><div class="">          Mat Object: (fieldsplit_1_) 1 MPI processes</div><div class="">            type: schurcomplement</div><div class="">            rows=147, cols=147</div><div class="">              Schur complement A11 - A10 inv(A00) A01</div><div class="">              A11</div><div class="">                Mat Object: 1 MPI processes</div><div class="">                  type: seqaij</div><div class="">                  rows=147, cols=147</div><div class="">                  total: nonzeros=147, allocated nonzeros=147</div><div class="">                  total number of mallocs used during MatSetValues calls=0</div><div class="">                    not using I-node routines</div><div class="">              A10</div><div class="">                Mat Object: 1 MPI processes</div><div class="">                  type: seqaij</div><div class="">                  rows=147, cols=512</div><div class="">                  total: nonzeros=2560, allocated nonzeros=2560</div><div class="">                  total number of mallocs used during MatSetValues calls=0</div><div class="">                    using I-node routines: found 87 nodes, limit used is 5</div><div class="">              KSP of A00</div><div class="">                KSP Object: (fieldsplit_0_) 1 MPI processes</div><div class="">                  type: gmres</div><div class="">                    restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div class="">                    happy breakdown tolerance 1e-30</div><div class="">                  maximum iterations=10000, initial guess is zero</div><div class="">                  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.</div><div class="">                  left preconditioning</div><div class="">                  using PRECONDITIONED norm type for convergence test</div><div class="">                PC Object: (fieldsplit_0_) 1 MPI processes</div><div class="">                  type: ilu</div><div class="">                    out-of-place factorization</div><div class="">                    0 levels of fill</div><div class="">                    tolerance for zero pivot 2.22045e-14</div><div class="">                    matrix ordering: natural</div><div class="">                    factor fill ratio given 1., needed 1.</div><div class="">                      Factored matrix follows:</div><div class="">                        Mat Object: 1 MPI processes</div><div class="">                          type: seqaij</div><div class="">                          rows=512, cols=512</div><div class="">                          package used to perform factorization: petsc</div><div class="">                          total: nonzeros=9213, allocated nonzeros=9213</div><div class="">                          total number of mallocs used during MatSetValues calls=0</div><div class="">                            not using I-node routines</div><div class="">                  linear system matrix = precond matrix:</div><div class="">                  Mat Object: 1 MPI processes</div><div class="">                    type: seqaij</div><div class="">                    rows=512, cols=512</div><div class="">                    total: nonzeros=9213, allocated nonzeros=9213</div><div class="">                    total number of mallocs used during MatSetValues calls=0</div><div class="">                      not using I-node routines</div><div class="">              A01</div><div class="">                Mat Object: 1 MPI processes</div><div class="">                  type: seqaij</div><div class="">                  rows=512, cols=147</div><div class="">                  total: nonzeros=2562, allocated nonzeros=2562</div><div class="">                  total number of mallocs used during MatSetValues calls=0</div><div class="">                    not using I-node routines</div><div class="">    linear system matrix = precond matrix:</div><div class="">    Mat Object: 1 MPI processes</div><div class="">      type: seqaij</div><div class="">      rows=659, cols=659</div><div class="">      total: nonzeros=14482, allocated nonzeros=27543</div><div class="">      total number of mallocs used during MatSetValues calls=1309</div><div class="">        not using I-node routines</div><br class=""><br class=""><br class=""><div class="gmail_quote_attribution">On Oct 9 2020, at 2:17 am, Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>> wrote:</div><blockquote class=""><div class=""><br class=""></div><div class="">  When you get a huge change at restart this means something is seriously wrong with either the linear operator or the linear preconditioner. </div><div class=""><br class=""></div><div class="">  How are you doing the matrix vector product?   Note both the operator and preconditioner must be linear operators for GMRES.</div><div class=""><br class=""></div><div class="">  FGMRES allows the preconditioner to be nonlinear. You can try</div><div class=""><br class=""></div><div class="">  -ksp_type fgmres -ksp_monitor_true_residual</div><div class=""><br class=""></div><div class="">   Barry</div><div class=""><br class=""><div class=""><br class=""><blockquote class=""><div class="">On Oct 8, 2020, at 2:43 AM, Yang Juntao <<a href="https://link.getmailspring.com/link/AA5B8AA2-FFAB-403A-A231-10CE902ACB44@getmailspring.com/0?redirect=mailto%3AY.Juntao%40hotmail.com&recipient=YnNtaXRoQHBldHNjLmRldg%3D%3D" title="mailto:Y.Juntao@hotmail.com" class="">Y.Juntao@hotmail.com</a>> wrote:</div><br class=""><div class=""><div class="WordSection1"><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">Hello, </font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class=""> </font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">I’m working on a nonlinear solver with SNES with handcoded jacobian and function. Each linear solver is solved with KSP solver.</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">But sometimes I got issues with ksp solver convergence. I tried with finite difference approximated jacobian, but get the same error.</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class=""> </font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">From the iterations, the convergence seems ok at the beginning but suddenly diverged in the last iteration.</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">Hope anyone with experience on ksp solvers could direct me to a direction I can debug the problem.</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class=""> </font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iter = 0, SNES Function norm 2.94934e-06</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 0 KSP Residual norm 1.094600281831e-06</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 1 KSP Residual norm 1.264284474186e-08</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 2 KSP Residual norm 6.593269221816e-09</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 3 KSP Residual norm 1.689570779457e-09</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 4 KSP Residual norm 1.040661505932e-09</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 5 KSP Residual norm 5.422761817348e-10</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 6 KSP Residual norm 2.492867371369e-10</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 7 KSP Residual norm 8.261522376775e-11</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 8 KSP Residual norm 4.246401544245e-11</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 9 KSP Residual norm 2.514366787388e-11</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 10 KSP Residual norm 1.982940267051e-11</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 11 KSP Residual norm 1.586470414676e-11</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 12 KSP Residual norm 9.866392216207e-12</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 13 KSP Residual norm 4.951342176999e-12</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 14 KSP Residual norm 2.418292660318e-12</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 15 KSP Residual norm 1.747418526086e-12</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 16 KSP Residual norm 1.094150535809e-12</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 17 KSP Residual norm 4.464287492066e-13</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 18 KSP Residual norm 3.530090494462e-13</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 19 KSP Residual norm 2.825698091454e-13</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 20 KSP Residual norm 1.950568425807e-13</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 21 KSP Residual norm 1.227898091813e-13</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 22 KSP Residual norm 5.411106347374e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 23 KSP Residual norm 4.511115848564e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 24 KSP Residual norm 4.063546606691e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 25 KSP Residual norm 3.677694771949e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 26 KSP Residual norm 3.459244943466e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 27 KSP Residual norm 3.263954971093e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 28 KSP Residual norm 3.087344619079e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 29 KSP Residual norm 2.809426925625e-14</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">iteration 30 KSP Residual norm 4.366149884754e-01</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  Linear solve did not converge due to DIVERGED_DTOL iterations 30</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class=""><strong class=""> </strong></font></font></font></font></div><div class=""><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class=""> </font></font></font></font></div></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">SNES Object: 1 MPI processes</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  type: newtonls</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  SNES has not been set up so information may be incomplete</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  maximum iterations=50, maximum function evaluations=10000</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  total number of linear solver iterations=0</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  total number of function evaluations=0</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  norm schedule ALWAYS</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  SNESLineSearch Object: 1 MPI processes</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    type: bt</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      interpolation: cubic</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      alpha=1.000000e-04</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    maxstep=1.000000e+08, minlambda=1.000000e-12</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    maximum iterations=40</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  KSP Object: 1 MPI processes</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    type: gmres</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      happy breakdown tolerance 1e-30</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    maximum iterations=10000, initial guess is zero</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    tolerances:  relative=1e-08, absolute=1e-50, divergence=10000.</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    left preconditioning</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    using DEFAULT norm type for convergence test</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">  PC Object: 1 MPI processes</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    type: fieldsplit</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    PC has not been set up so information may be incomplete</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      FieldSplit with Schur preconditioner, factorization FULL</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      Preconditioner for the Schur complement formed from S itself</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      Split info:</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      KSP solver for A00 block</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">          not yet available</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      KSP solver for S = A11 - A10 inv(A00) A01</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">          not yet available</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    linear system matrix = precond matrix:</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">    Mat Object: 1 MPI processes</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      type: seqaij</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      rows=659, cols=659</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      total: nonzeros=659, allocated nonzeros=7908</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">      total number of mallocs used during MatSetValues calls=0</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">        not using I-node routines</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class=""> </font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">Regards</font></font></font></font></div><div class=""><font style="font-size:18px" class=""><font style="font-family:Helvetica" class=""><font style="font-size:11pt" class=""><font style="font-family:Calibri, sans-serif" class="">Juntao</font></font></font></font></div></div></div></blockquote></div></div></blockquote><img class="mailspring-open" alt="Sent from Mailspring" width="0" height="0" style="border:0; width:0; height:0;" src="https://link.getmailspring.com/open/AA5B8AA2-FFAB-403A-A231-10CE902ACB44@getmailspring.com?me=05842e25&recipient=YnNtaXRoQHBldHNjLmRldg%3D%3D"></div></blockquote></div><br class=""></body></html>