<div dir="ltr"><div dir="ltr">On Fri, May 27, 2022 at 9:42 AM Sami BEN ELHAJ SALAH <<a href="mailto:sami.ben-elhaj-salah@ensma.fr">sami.ben-elhaj-salah@ensma.fr</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 style="overflow-wrap: break-word;"><div><div>Hello Isaac,</div><div><br></div><div>Thank you for your reply!</div><div><br></div><div>Let me confirm that when I use DMCreateMatrix() with the orig_dm, I got my jacobian_matrix. Also, I have succeeded to solve my system and my solution was converged<b> by using one process</b>.</div><div>Let me give you some other information about my code. Currently, I am using my own discretization system and not the PetscDS object. Considering the nonlinear solver SNES, I use the following routines (like the basic snes usage):</div><div>- SNESSetFunction(snes,residual_vector,compute_residual ,(void*) _finite_element_formulation)</div><div>- SNESSetJacobian(snes, jacobian.matrix(), _jacobian_precond_matrix, compute_jacobian, (void*) _finite_element_formulation) </div><div>- SNESSolve(snes,NULL,x)</div><div><br></div><div>Regarding  your last answer, I will try to reformulate my question as follows:</div><div>Using a distributed dm instead of the original dm (not distributed dm) and my own discretization system (not the PetscDS object),</div></div></div></blockquote><div><br></div><div>You do not have to use PetscDS, but the DM does need to a PetscSection in order to compute sizes and sparsity patterns. Do you have that?</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 style="overflow-wrap: break-word;"><div><div> should I add something specific to get a distributed jacobian_matrix over processes?</div><div>I precise that I just replaced the orig_dm by a distributed mesh with the routine that I wrote in my first mail. So is it enough ?</div><div> </div><div>Thank you and have a good day,</div><div>Sami,</div></div><div><br></div><br><div>
<div dir="auto" style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div dir="auto" style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div dir="auto" style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div dir="auto" style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none">--</div><div dir="auto" style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none">Dr. Sami BEN ELHAJ SALAH<br>Ingénieur de Recherche (CNRS)<br>Institut Pprime - ISAE - ENSMA<br>Mobile: 06.62.51.26.74<br><a href="mailto:sami.ben-elhaj-salah@ensma.fr" target="_blank">Email: sami.ben-elhaj-salah@ensma.fr</a></div><div dir="auto" style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><a href="https://samiben91.github.io/samibenelhajsalah/index.html" target="_blank">www.samibenelhajsalah.com</a><br><br><br></div></div></div></div>
</div>
<div><br><blockquote type="cite"><div>Le 25 mai 2022 à 19:41, Toby Isaac <<a href="mailto:toby.isaac@gmail.com" target="_blank">toby.isaac@gmail.com</a>> a écrit :</div><br><div><div>Hi Sami,<br><br>Just to verify: if you call DMCreateMatrix() on orig_dm, do you get a<br>Jacobian matrix?<br><br>The DMPlex must be told what kind of discretized fields you want a<br>matrix for and what equations you are discretizing.  This is handled<br>by the PetscDS object.  In snes/tutorials/ex59.c, see the code after<br>DMGetDS() for an example.<br><br>- Toby<br><br>On Wed, May 25, 2022 at 1:17 PM Sami BEN ELHAJ SALAH<br><<a href="mailto:sami.ben-elhaj-salah@ensma.fr" target="_blank">sami.ben-elhaj-salah@ensma.fr</a>> wrote:<br><blockquote type="cite"><br>Dear PETSc developer team,<br><br>I m trying to create à jacobian_matrix from distributed DM. I have followed the two examples (snes/tests/ex2.c and ex56.c). So I wrote this routine:<br><br>PetscDM orig_dm;<br>PetscDM dist_dm = NULL;<br>PetscPartitioner part;<br>DMPlexGetPartitioner(orig_dm, &part);<br>PetscPartitionerSetType(part, PETSCPARTITIONERPARMETIS);<br>DMPlexDistribute(orig_dm, 0, NULL, &dist_dm);<br><br>PetscErrorCode err = DMCreateMatrix(dist_dm, &jacobian_matrix);<br>PetscInt M, N, m, n;<br>MatGetSize(jacobian_matrix, &M, &N);<br>MatGetLocalSize(jacobian_matrix, &m, &n);<br><br>Then I run my code with 2 processes and I obtained this result:<br>Size from jacobian_matrix: M=0 m =0 (this result is the same in all processes).<br><br>I didn't understand if I forgot something in my code to obtain the correct values for the local and global sizes for the jacobean_matrix? (I mean if my code is missing something to obtain a distributed mesh over processes ?)<br><br>Thank you in advance for any help!<br>Sami,<br><br><br>--<br>Dr. Sami BEN ELHAJ SALAH<br>Ingénieur de Recherche (CNRS)<br>Institut Pprime - ISAE - ENSMA<br>Mobile: 06.62.51.26.74<br><a href="mailto:sami.ben-elhaj-salah@ensma.fr" target="_blank">Email: sami.ben-elhaj-salah@ensma.fr</a><br><a href="http://www.samibenelhajsalah.com" target="_blank">www.samibenelhajsalah.com</a><br><br><br><br></blockquote></div></div></blockquote></div><br></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>