<div dir="ltr"><div><div><div>Thanks Barry and Matt,<br><br></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><br></div>Cheers,<br></div>Saransh<br></div><br><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">bsmith@petsc.dev</a>> wrote:<br></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;"><div> Saransh,</div><div> </div> I've add some code for SNESSetPicard() in the PETSc branch barry/2021-05-06/add-snes-picard-mf see also http<a>s://gitlab.com/petsc/petsc/-/merge_requests/3962</a> that will make your coding much easier.<div><br></div><div> With this branch you can provide code that computes A(x), using SNESSetPicard(). </div><div><br></div><div>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><br></div><div>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></div><div><br></div><div>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><br></div><div>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><br></div><div>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><br></div><div> Hope this helps,</div><div><br></div><div> Barry</div><div><br></div><div><div><div><br><blockquote type="cite"><div>On May 6, 2021, at 1:21 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div dir="ltr">On Thu, May 6, 2021 at 2:09 PM Saransh Saxena <<a href="mailto:saransh.saxena5571@gmail.com" target="_blank">saransh.saxena5571@gmail.com</a>> wrote:<br></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">Hi,<br><br>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><br>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><br></div><div>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></div></div></blockquote><div><br></div><div>You would give PETSc an outer function MyResidual() that looked like this:</div><div><br></div><div>PetscErrorCode MyResidual(SNES snes, Vec X, Vec F, void *ctx)</div><div>{</div><div> <call your code to compute b, or pass it in using ctx></div><div> <call your code to compute A(X)></div><div> MatMult(A, X, F);</div><div> VecAXPY(F, -1.0, b);</div><div>}</div><div> </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"><div>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></div></div></blockquote><div><br></div><div>If you do nothing, we will compute it by default.</div><div><br></div><div> Thanks,</div><div><br></div><div> MAtt</div><div> </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"><div>Best regards,<br><br>Saransh</div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</div></blockquote></div><br></div></div></div></blockquote></div>