<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=""><blockquote type="cite" class=""><div dir="ltr" class=""><div class=""><b class="">ierr = MatSeqAIJGetArray(ptr->K,&Kc);CHKERRQ(ierr); </b></div><div class=""><b class="">ierr = MatGetLocalSize(ptr->K,&m,&n); CHKERRQ(ierr);<br class=""></b></div><b class="">for(i=0;i<m;i++)<br class="">    { <br class="">      for(j=0;j<n;j++)<br class="">      {<br class="">        v  = *(Kc + (i*n+j));<br class="">        ierr = MatSetValues(Kloc,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);<br class="">      }<br class=""></b><div class=""><b class="">    } </b><br class=""><div class="">    </div></div></div></blockquote><div class=""><br class=""></div>I don't understand the purpose of this code fragment. Are you wanting to copy your K into Kloc? What is the purpose of Kloc? With standard usage of SNES one only needs to provide PETSc a global matrix, one does not need to work with "ghosted" matrices at all. <div class=""><br class=""></div><div class="">Also normally one fills the f vector in FormFunction and the Jacobian in FormJacobian and there is no need to access the Jacobians at all in SNES <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html#functionjacobian" class="">https://www.mcs.anl.gov/petsc/documentation/faq.html#functionjacobian</a>  nor do you need to pass into FormFunction in the context the F or Jacobian matrix. </div><div class=""><br class=""></div><div class="">Maybe if you explained your use case a bit more we could make suggestions on how to accomplish it.</div><div class=""><br class=""></div><div class="">Barry</div><div class=""><br class=""></div><div class=""><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On May 19, 2021, at 10:00 AM, Kaushik Vijaykumar <<a href="mailto:kaushikv318@gmail.com" class="">kaushikv318@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class="">Hi everyone,<div class=""><br class=""></div><div class="">We are in the process of integrating nonlinear solution to our FE code through SNES. The first steps that I understood that need to be done is to be able to pass the assembled stiffness matrix K and force vector F to the "<b class="">formfunction</b>" to calculate the residual and "<b class="">formjacobian</b>" to calculate the tangent stiffness matrix. To do so, I have defined a struct:</div><div class=""><br class=""></div><div class=""><b class="">extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);<br class="">extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);</b><br class=""><br class="">// define a struct type to pass K and f as "user context"<br class=""><b class="">typedef struct<br class="">    {<br class="">      Mat K;<br class="">      Vec f;<br class="">    } K_and_f;</b><br class=""></div><div class=""><br class=""></div><div class="">In the main program, the struct is declared</div><div class=""><b class="">int main()</b></div><div class=""><b class="">{</b></div><div class=""><b class="">K_and_f main_K_and_f;</b> // declare the struct</div><div class=""><br class=""></div>// SNES - Populate K and f into the struct<br class=""><b class="">main_K_and_f.K = K;</b> // K matrix<br class=""><div class=""><b class="">main_K_and_f.f = f; </b>// f vector <br class=""></div><div class="">....</div><div class="">....</div><div class=""><b class="">}</b></div><div class="">In form function</div><div class=""><br class=""></div><div class=""><b class="">PetscErrorCode FormFunction(SNES snes, Vec x, Vec F, void* ctx) {<br class="">    PetscErrorCode ierr;<br class="">    PetscReal *ax,*c;<br class="">    PetscReal *aF;<br class="">    PetscScalar *Kc,v; <br class="">    PetscInt nlocal,m,n,i,j,index[100000];<br class=""><br class="">   </b> // Create local Vec and Mat<br class=""><b class="">    Mat Kloc;</b><br class=""><b class="">    Vec Floc, Uloc, KUloc, resloc;</b><br class=""><b class="">    </b><br class=""><b class="">    K_and_f* ptr = (K_and_f*) ctx; </b>// cast the pointer to void into pointer to struct<br class=""></div><div class=""><b class="">    </b>// Get local F array, FLOC</div><b class="">    ierr = VecGhostGetLocalForm(ptr->f,&Floc);CHKERRQ(ierr);</b><div class=""><br class=""></div><div class=""><br class=""></div><div class="">I am able to get the f array from the main program using "

VecGhostGetLocalForm" and the vec is correct. However, I am having trouble with reading the K matrix in formfunction. </div><div class=""><br class=""></div><div class="">The following trial gives me incorrect K values in Kloc:</div><div class=""><b class="">ierr = MatSeqAIJGetArray(ptr->K,&Kc);CHKERRQ(ierr); </b></div><div class=""><b class="">ierr = MatGetLocalSize(ptr->K,&m,&n); CHKERRQ(ierr);<br class=""></b></div><b class="">for(i=0;i<m;i++)<br class="">    { <br class="">      for(j=0;j<n;j++)<br class="">      {<br class="">        v  = *(Kc + (i*n+j));<br class="">        ierr = MatSetValues(Kloc,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);<br class="">      }<br class=""></b><div class=""><b class="">    } </b><br class=""><div class="">    </div><div class="">    When I compare K and Kloc, they are not identical </div><div class=""><br class=""></div><div class=""><br class=""></div><div class="">Please let me know if there is an equivalent function like <b class="">VecGhostGetLocalForm()</b> for Matrices. If not, is there a better way to do this.</div><div class=""><br class=""></div><div class="">Any guidance/help is greatly appreciated.</div><div class=""><br class=""></div><div class="">Thanks</div><div class="">Kaushik</div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div></div></div>
</div></blockquote></div><br class=""></div></div></body></html>