<html><head><meta http-equiv="Content-Type" content="text/html; charset=us-ascii"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div>  Saransh,<div class=""><br class=""></div><div class="">     If Picard or Newton's method does not converge, you can consider adding pseudo-transient and/or other continuation methods. For example, if the problem is made difficult by certain physical parameters you can start with "easier" values of the parameters, solve the nonlinear system, then use its solution as the initial guess for slightly more "difficult" parameters, etc. Or, depending on the problem grid sequencing may be appropriate. We have some tools to help with all these approaches.</div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On May 9, 2021, at 2:07 PM, Saransh Saxena <<a href="mailto:saransh.saxena5571@gmail.com" class="">saransh.saxena5571@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class=""><div class=""><div class=""><div class="">Thanks Barry and Matt,<br class=""><br class=""></div>Till now I was only using a simple fixed point nonlinear solver manually coded instead of ones provided by PETSc. However, the problem I am trying to solve is highly nonlinear so I suppose I'll need at least a newton based solver to start with. I'll get back to you guys if I have any questions.<br class=""><br class=""></div>Cheers,<br class=""></div>Saransh<br class=""></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sat, May 8, 2021 at 5:18 AM Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>> wrote:<br 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 style="overflow-wrap: break-word;" class=""><div class="">  Saransh,</div><div class=""> </div>  I've add some code for SNESSetPicard() in the PETSc branch barry/2021-05-06/add-snes-picard-mf see also http<a class="">s://gitlab.com/petsc/petsc/-/merge_requests/3962</a> that will make your coding much easier.<div class=""><br class=""></div><div class="">   With this branch you can provide code that computes A(x), using SNESSetPicard(). </div><div class=""><br class=""></div><div class="">1) by default it will use the defection-correction form of Picard iteration  A(x^n)(x^{n+1} - x^{n}) = b - A(x^n)  to solve, which can be slower than Newton  </div><div class=""><br class=""></div><div class="">2) with -snes_fd_color it will compute the Jacobian via coloring using SNESComputeJacobianDefaultColor() (assuming the true Jacobian has the same sparsity structure as A). The true Jacobian is J(x^n) = A'(x^n)[x^n] - A(x^n) where A'(x^n) is the third order tensor of the derivatives of A() and A'(x^n)[x^n] is a matrix, I do not know if, in general, it has the same nonzero structure as A. (I'm lost beyond matrices :-().</div><div class=""></div><div class=""><br class=""></div><div class="">3) with -snes_mf_operator it will apply the true Jacobian matrix-free and precondition it with a preconditioner built from A(x^n) matrix, for some problems this works well. </div><div class=""><br class=""></div><div class="">4) with -snes_fd it uses SNESComputeJacobianDefault() and computes the Jacobian by finite differencing one column at a time, thus it is very slow and not useful for large problems. But useful for testing with small problems.</div><div class=""><br class=""></div><div class="">So you can provide A() and need not worrying about providing the Jacobian or even the function evaluation code. It is all taken care of by SNESSetPicard().</div><div class=""><br class=""></div><div class="">  Hope this helps,</div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""></div><div class=""><div class=""><div class=""><br class=""><blockquote type="cite" class=""><div class="">On May 6, 2021, at 1:21 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank" class="">knepley@gmail.com</a>> wrote:</div><br class=""><div class=""><div dir="ltr" class=""><div dir="ltr" class="">On Thu, May 6, 2021 at 2:09 PM Saransh Saxena <<a href="mailto:saransh.saxena5571@gmail.com" target="_blank" class="">saransh.saxena5571@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="">Hi,<br class=""><br class="">I am trying to incorporate newton method in solving a nonlinear FEM equation using SNES from PETSc. The overall equation is of the type A(x).x = b, where b is a vector of external loads, x is the solution field (say displacements for e.g.) and A is the combined LHS matrix derived from the discretization of weak formulation of the governing finite element equation. <br class=""><br class="">While going through the manual and examples of snes, I found that I need to define the function of residual using SNESSetFunction() and jacobian using SNESSetJacobian(). In that context I had a couple of questions :-<div class=""><br class=""></div><div class="">1. In the snes tutorials I've browsed through, the functions for computing residual passed had arguments only for x, the solution vector and f, the residual vector. Is there a way a user can pass an additional vector (b) and matrix (A) for computing the residual as well? as in my case, f = b - A(x).x<br class=""></div></div></blockquote><div class=""><br class=""></div><div class="">You would give PETSc an outer function MyResidual() that looked like this:</div><div class=""><br class=""></div><div class="">PetscErrorCode MyResidual(SNES snes, Vec X, Vec F, void *ctx)</div><div class="">{</div><div class="">  <call your code to compute b, or pass it in using ctx></div><div class="">  <call your code to compute A(X)></div><div class="">  MatMult(A, X, F);</div><div class="">  VecAXPY(F, -1.0, b);</div><div class="">}</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=""><div class="">2. Since computing jacobian is not that trivial, I would like to use one of the pre-built jacobian methods. Is there any other step other than setting the 3rd argument in SNESSetJacobian to SNESComputeJacobianDefault?<br class=""></div></div></blockquote><div class=""><br class=""></div><div class="">If you do nothing, we will compute it by default.</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=""><div class="">Best regards,<br class=""><br class="">Saransh</div></div>
</blockquote></div><br clear="all" class=""><div class=""><br class=""></div>-- <br class=""><div dir="ltr" class=""><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></blockquote></div>
</div></blockquote></div><br class=""></div></body></html>