<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=""><br class=""></div><div class="">  Pranay</div><div class=""><br class=""></div>   Newton's method is generally the best choice for nonlinear problems as Matt notes but PETSc also provides an implementation of Picard's method with SNESSetPicard() <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetPicard.html" class="">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetPicard.html</a><div class=""><br class=""></div><div class="">   We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then the direct Picard iteration A(x^n) x^{n+1} = b(x^n), which is what Matt just said.</div><div class=""><br class=""></div><div class="">   Based on your email it looks like you using the direct Picard iteration algorithm. </div><div class=""><br class=""></div><div class="">   With your current code you can likely easily switch to trying SNESSetPicard() and then switch to trying Newton with  SNESSetFunction(), <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html" class="">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html</a>  and SNESSetJacobian() <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html" class="">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html</a></div><div class=""><br class=""></div><div class="">   PETSc is designed to make it as simple as possible to switch between various algorithms to help determine the best for your exact problem. </div><div class=""><br class=""></div><div class="">   SNES uses KSP for is linear solvers so you get full access to all the possible preconditioners, for larger problems as Matt notes, once you have the best nonlinear convergence selected tuning the linear solver is the most important thing to do for speed. We recommend when possible first getting good nonlinear convergence using a direct linear solver and then switching to an iterative solver as an optimization, for large problems you almost always should to switch to an iterative solver when the problem size increases.<br class=""><div class=""><br class=""></div><div class="">    Barry</div><div class=""><br class=""><div class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Aug 22, 2020, at 7:24 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class=""><div dir="ltr" class="">On Sat, Aug 22, 2020 at 2:07 AM baikadi pranay <<a href="mailto:pranayreddy865@gmail.com" class="">pranayreddy865@gmail.com</a>> wrote:<br class=""></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr" class="">Hello,<br class=""><br class="">I am trying to solve the Poisson equation in 2D for heterostructure devices. I have linearized the equation and discretized it using FDM. I am using BiCGStab to iteratively solve for the solution as follows:<br class=""><br class="">Step 1: Solve A^(i-1) x^(i) = b^(i-1)     {i = 1 to N where convergence is reached}<br class="">Step 2: Use x^{i} to update the central coefficients of A^{i-1} to get A^{i} and                          similarly  update b^{i-1} to get b^{i}<br class="">Step3: If ( ||x^{i}-x^{i-1}||_2 , the 2-norm of the solution update, is greater than a                     tolerance, then go back to Step 1  to solve the new system of equations                     using BiCGStab. Else, exit the loop.<br class=""><b class=""><u class="">1) I am facing the following problem with this procedure</u></b>:<br class="">The 2-norm of the solution update is suddenly becoming zero after a few iterations in some cases. I print out the getconvergedreason and there are not red flags there, so I am kind of confused whey this behaviour is being observed. This behaviour is leading to "false convergences", in the sense that the solutions obtained are not physical.<br class=""><br class="">A similar behaviour was observed when I used SOR instead of BiCGStab. At this point I am starting to suspect if it is wrong to use linear solvers on the poisson equation which is a nonlinear equation (although linearized). If you could please comment on this, that would be very helpful.  <br class=""><br class="">Any help with this problem is greatly appreciated. Please let me know if you need any further information. <br class=""></div></blockquote><div class=""><br class=""></div><div class="">1) You are coding up the Picard method by hand to solve your nonlinear equation. If the operator is not contractive, this can stagnate, as you are seeing. You</div><div class="">    could try another solver, like Newton's method. We have a variety of nonlinear solves in the SNES class.</div><div class=""><br class=""></div><div class="">2) It is not clear from your description whether you linear solver is converging. BiCGStab without a preconditioner is a terrible solver for Poisson. We usually</div><div class="">    recommend starting with Algebraic Multigrid, like Hypre which is great at 2D Poisson. You can monitor the convergence of your linear solver using</div><div class=""><br class=""></div><div class="">      -knp_monitor_true_solution -ksp_converged_reason</div><div class=""><br class=""></div><div class="">   We want to see this information with any questions about convergence.</div><div class=""><br class=""></div><div class="">  Thanks,</div><div class=""><br class=""></div><div class="">     Matt</div><div class=""> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr" class="">Thank you,<br class=""><br class="">Sincerely,<br class="">Pranay.<br class=""><br class=""></div><div hspace="streak-pt-mark" style="max-height:1px" class=""><img alt="" style="width: 0px; max-height: 0px; overflow: hidden;" src="https://mailfoogae.appspot.com/t?sender=acHJhbmF5cmVkZHk4NjVAZ21haWwuY29t&type=zerocontent&guid=e91af943-626d-449a-b812-eafec67b8f5a" class=""><font color="#ffffff" size="1" class="">ᐧ</font></div>
</blockquote></div><br clear="all" class=""><div class=""><br class=""></div>-- <br class=""><div dir="ltr" class="gmail_signature"><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">-- Norbert Wiener</div><div class=""><br class=""></div><div class=""><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank" class="">https://www.cse.buffalo.edu/~knepley/</a><br class=""></div></div></div></div></div></div></div></div>
</div></blockquote></div><br class=""></div></div></div></body></html>