<div dir="ltr"><div dir="ltr">On Wed, Feb 16, 2022 at 8:18 AM Bojan Niceno <<a href="mailto:bojan.niceno.scientist@gmail.com">bojan.niceno.scientist@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">Dear Matthew and Barry,<div><br></div><div>I am going through tutorial examples, and have noticed that calls to <font face="monospace">MatMPIAIJSetPreallocation</font> and <font face="monospace">MatSeqAIJSetPreallocation</font> always come in succession of one another, usually in that order.  Why would that be?  Is it maybe to be able to run on systems without MPI, so that the first function <font face="monospace">MatMPIAIJSetPreallocation</font> only makes a dummy call?</div><div><br></div><div>I was checking the manual pages and found no information about this pairing.</div></div></blockquote><div><br></div><div>This is because we are too lazy or time-constrained to convert all examples to the new style using</div><div><br></div><div>  <a href="https://petsc.org/main/docs/manualpages/Mat/MatXAIJSetPreallocation.html">https://petsc.org/main/docs/manualpages/Mat/MatXAIJSetPreallocation.html</a></div><div><br></div><div>More broadly, PETSc is constructed like Objective-C in that method bindings can change at runtime and calling a</div><div>method on an object of a different type just ignores it. So if I have a SeqAIJ matrix, only the <span style="font-family:monospace">MatSeqAIJSetPreallocation</span></div><div>call works and the other is ignored. This way you can put in customization for many types without lots of checks</div><div>to see what type the caller has.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div><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 dir="ltr"><div><br></div><div>    Kind regards,</div><div><br></div><div>    Bojan</div><div><br></div><div><br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Feb 16, 2022 at 7:30 AM Bojan Niceno <<a href="mailto:bojan.niceno.scientist@gmail.com" target="_blank">bojan.niceno.scientist@gmail.com</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 dir="ltr">Thanks you guys for your suggestions, I will have a look at functions you suggested.<div><br></div><div>    Have a great day,</div><div><br></div><div>    Bojan</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Tue, Feb 15, 2022 at 11:29 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"><div><div><br></div>MatCreateMPIAIJWithArrays() and MatUpdateMPIAIJWithArrays() may be suitable for your use case.<div><br></div><div>This should also be much more efficient in moving the matrix from your code to PETSc's format.<br><div><br></div><div>  Barry</div><div><br><div><br><blockquote type="cite"><div>On Feb 15, 2022, at 4:13 PM, Bojan Niceno <<a href="mailto:bojan.niceno.scientist@gmail.com" target="_blank">bojan.niceno.scientist@gmail.com</a>> wrote:</div><br><div><div dir="ltr">Dear PETSc users,<div><br></div><div><br></div><div>I have an in-house computational fluid dynamics (CFD) solver, written in Fortran 2008, parallelized with MPI with its own home-grown suite of linear solvers.  The code is unstructured, performs domain decomposition with METIS and all communication buffers, I mean connectivity between processors, has been properly worked out.</div><div><br></div><div>A couple of weeks back, I decided to try out the PETSc suite of solvers.  After some initial setbacks, I managed to compile my code with PETSc and have the sequential version working fine :-)</div><div><br></div><div>I have essentially using the following PETSc routines to get the code solving linear systems with PETSc:  </div><div><br></div><div><i> </i>I set up the working space as follows:</div><div><br></div><div><font face="monospace">  call PetscInitialize(PETSC_NULL_CHARACTER, Pet % petsc_err)<br></font></div><div><font face="monospace">  call MatCreateSeqAij(PETSC_COMM_SELF, ...<br></font></div><div><font face="monospace">  </font><font face="monospace" color="#ff0000">call MatSeqAijSetColumnIndices(....</font></div><div><font face="monospace">  call VecCreateSeq(PETSC_COMM_SELF,  ...  ! Create Vector x</font></div><div><font face="monospace">  call VecCreateSeq(PETSC_COMM_SELF,  ...  ! Create Vector b</font></div><div><font face="monospace">  call KSPCreate(PETSC_COMM_SELF,  ...</font></div><div><font face="monospace"><br></font></div><div><font face="arial, sans-serif">Then in order to solve a system, I do:</font></div><div><font face="arial, sans-serif"><br></font></div><div><font face="monospace">  call MatSetValue(Pet % petsc_A,   ! Inside a loop through matrix entries<br>                   i-PETSC_ONE, <br>                   k-PETSC_ONE, ...</font></div><font face="monospace">  call MatAssemblyBegin(Pet % petsc_A, MAT_FINAL_ASSEMBLY, Pet % petsc_err)<br>  call MatAssemblyEnd  (Pet % petsc_A, MAT_FINAL_ASSEMBLY, Pet % petsc_err)<br></font><div><font face="monospace"> <br></font></div><div><font face="monospace">  call VecSetValue(Pet % petsc_x, ... ! Fill up x</font></div><div><div><font face="monospace">  call VecSetValue(Pet % petsc_b, ... ! Fill up b</font></div><font face="monospace"><div><font face="monospace"><br></font></div>  call KSPSetType(Pet % petsc_ksp ... ! Set solver<br></font></div><div><font face="monospace">  call KSPGetPC(Pet % petsc_ksp, ...  ! get preconditioner context</font></div><font face="monospace">  call PCSetType(Pet % petsc_pc, ...  ! Set preconditioner</font><div><font face="monospace"><br>  call KSPSetFromOptions(Pet % petsc_ksp, Pet % petsc_err)                       <br>  call KSPSetUp         (Pet % petsc_ksp, Pet % petsc_err)                       <br></font><div><font face="monospace"> </font></div><div><font face="monospace">  ! Finally solve</font></div><div><font face="monospace">  call KSPSolve(Pet % petsc_ksp, ...</font></div><div><font face="monospace"><br></font></div><div><font face="arial, sans-serif">Once this was up and running, I thought that in order to have the parallel version I will merely have to replace the "Seq" versions of the above functions, with their parallel counterparts.  I was expecting to find the red function (</font><font face="monospace">MatSeqAijSetColumnIndices</font><font face="arial, sans-serif">) for parallel runs, but it doesn't seem to exist.  I have found non-seq versions of some other functions (</font><font face="monospace">MatCreateAij</font><font face="arial, sans-serif">, </font><font face="monospace">VecCreateSeq</font><font face="arial, sans-serif">), but not something like </font><font face="monospace">MatAijSetColumnIndices</font><font face="arial, sans-serif">, which surprised me a bit, because I have this information in my code.</font></div><div><span style="font-family:arial,sans-serif">  </span></div><div><div><font face="arial, sans-serif">Is there a parallel counterpart of this function, and if there is none, what should it be replaced with?  I understand that I will have to provide non-zeros in buffers (o_nnz), which is not a big issue, but how to provide information on columns for parallel version is not clear to me.  In a nutshell, I would need a hint on which of the above functions could remain the same in parallel, and which should be replaced and with what?</font></div><div><br></div></div></div><div>    Cheers,</div><div><br></div><div>    Bojan</div><div><br></div></div>
</div></blockquote></div><br></div></div></div></blockquote></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>