<div dir="ltr"><div dir="ltr">On Mon, Jun 14, 2021 at 7:21 AM Kaushik Vijaykumar <<a href="mailto:kaushikv318@gmail.com">kaushikv318@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">Thanks for the suggestion Barry, I did do that by setting it in multiple ways first by setting<div><br></div><div><div>Formfunction(snes,Vec x,Vec F,void* ctx)</div><div>{</div><div>ierr = SNESGetJacobian(snes,&jac,NULL,NULL,NULL);CHKERRQ(ierr);<br></div><div>ierr = MatZeroEntries(jac);CHKERRQ(ierr);<br></div></div><div><br></div><div>.....</div><div><br></div><div><br></div><div>ierr = MatGetOwnershipRange(jac,&istart,&iend);<br><br> for(i=istart;i<iend;i++)<br> {<br> val = 0.0;<br> ierr = MatSetValue(jac, i, i, val, ADD_VALUES);CHKERRQ(ierr);<br> }<br></div><div>ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br> ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br></div><div>}</div><div>I get the same error as before for a single processor run.</div><div>"[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------</div>[0]PETSC ERROR: Object is in wrong state<br>[0]PETSC ERROR: Matrix is missing diagonal entry 0"</div></blockquote><div><br></div><div>It must be that the matrix you are setting the diagonal for is not the same matrix as referenced in the error message.</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>Thanks</div><div>Kaushik</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Jun 11, 2021 at 11:33 AM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">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><br><div><blockquote type="cite"><div dir="ltr"><div>When I execute this I get a Petsc error, object in wrong state missing diagonal entry.</div></div></blockquote><div><br></div> When do you get this error? During a TSSolve? While KSP is building a preconditioner.</div><div><br></div><div> Some parts of PETSc require that you put entries (even if they are zero) on all diagonal entries in the Jacobian. So likely you need to just make sure that your populate jac routine puts/adds 0 to all diagonal entries.</div><div><br></div><div> Barry</div><div><br></div><div><br></div><div><br><blockquote type="cite"><div>On Jun 11, 2021, at 10:16 AM, Kaushik Vijaykumar <<a href="mailto:kaushikv318@gmail.com" target="_blank">kaushikv318@gmail.com</a>> wrote:</div><br><div><div dir="ltr">Thanks Barry for your comments. Based on your suggestion I tried the following<div><br></div><div>Main()</div><div>{</div><div>ierr = SNESSetFunction(snes,r1,FormFunction,&this_ctx);CHKERRQ(ierr);</div><div>ierr = SNESSetJacobian(snes,jac,NULL,FormJacobian,&this_ctx);CHKERRQ(ierr);<br> ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);<br> ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);<br></div><div>}</div><div><br></div><div>In formfunction I populate the jacobian "<b>jac</b>"</div><div><br></div><div>Formfunction(snes,Vec x,Vec F,void* ctx)</div><div>{</div><div>ierr = SNESGetJacobian(snes,&jac,Null,NULL,NULL);CHKERRQ(ierr);<br></div><div><br></div><div>// "<b>Populate jac</b>"</div><div>ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br> ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); <br></div><div><br></div><div>}</div><div><div>FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)</div><div>{</div><div>ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br> ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br><br> ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br> ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br></div><div>}</div></div><div><br></div><div>When I execute this I get a Petsc error, object in wrong state missing diagonal entry.</div><div><br></div><div>I am not sure what I am missing here? </div><div><br></div><div>Thanks</div><div>Kaushik</div><div><br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Thu, Jun 10, 2021 at 6:37 PM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">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"><br>
<br>
> On Jun 10, 2021, at 3:29 PM, Kaushik Vijaykumar <<a href="mailto:kaushikv318@gmail.com" target="_blank">kaushikv318@gmail.com</a>> wrote:<br>
> <br>
> Hi everyone,<br>
> <br>
> I am trying to copy the stiffness matrix that I generated in form function to the jacobian matrix jac in form jacobian. The piece of code for that:<br>
> <br>
> <br>
If you have created jac and B in your main program and passed them with SNESSetJacobian() then you should be able to just use<br>
<br>
> ierr = MatCopy(ptr->K,jac,SAME_NONZERO_PATTERN);CHKERRQ(ierr);<br>
<br>
Generally jac and B are the same matrix so you don't need the second copy. <br>
> ierr = MatCopy(ptr->K,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); <br>
<br>
If they are sometimes the same and sometime not you can do<br>
<br>
if (jac != B) {ierr = MatCopy(ptr->K,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); }<br>
<br>
<br>
The following won't work because you are creating new local matrices in your form jacobian that don't exist outside of that routine.<br>
<br>
> <br>
> ierr = MatDuplicate(ptr->K,MAT_COPY_VALUES,&jac);CHKERRQ(ierr);<br>
> ierr = MatDuplicate(ptr->K,MAT_COPY_VALUES,&B);CHKERRQ(ierr);<br>
> <br>
> ierr = MatCopy(ptr->K,jac,SAME_NONZERO_PATTERN);CHKERRQ(ierr);<br>
> ierr = MatCopy(ptr->K,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); <br>
> <br>
> <br>
If the matrix you are creating in form function is the Jacobian then you don't need to use the memory and time to make copies. You can just form the matrix directly into the Jacobian you will use for SNES. You can obtain this with a call <br>
<br>
SNESGetJacobian(snes,&jac,NULL,NULL,NULL); <br>
<br>
in your form function and just put the matrix values directly into jac. Then your form Jacobian does not need to do anything but return since the jac already contains the Jacobian.<br>
<br>
On the off-change you are using the iteration A(x^n) x^{n+1} = b(x^n} and not using Newtons' method (which is why you compute A in form function, then you might be better off using SNESSetPicard() instead of SNESSetFunction(), SNESSetJacobian().<br>
<br>
> <br>
> Thanks<br>
> Kaushik<br>
<br>
</blockquote></div>
</div></blockquote></div><br></div></blockquote></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><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>