<html>
  <head>
    <meta http-equiv="Content-Type" content="text/html; charset=utf-8">
  </head>
  <body text="#000000" bgcolor="#FFFFFF">
    <p>Hi, thanks for the fix. It seems to work fine.</p>
    <p>Out of curiosity, I noticed the MatOrderingType of
      DMPlexGetOrdering is not used. Is this intentional ? To match
      MatGetOrdering ?</p>
    <p><br>
    </p>
    <p>Pierre<br>
    </p>
    <br>
    <div class="moz-cite-prefix">On 27/10/21 15:03, Matthew Knepley
      wrote:<br>
    </div>
    <blockquote type="cite"
cite="mid:CAMYG4Gk7jMi8CakkWtEbYKqJ2G3dPX=GpBSDGtja-M3UUWVphQ@mail.gmail.com">
      <div dir="ltr">
        <div dir="ltr">On Wed, Oct 27, 2021 at 3:15 AM Pierre Seize <<a
            href="mailto:pierre.seize@onera.fr" moz-do-not-send="true">pierre.seize@onera.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 bgcolor="#FFFFFF">
              <p><br>
              </p>
              <br>
              <div>On 26/10/21 22:28, Matthew Knepley wrote:<br>
              </div>
              <blockquote type="cite">
                <div dir="ltr">
                  <div dir="ltr">On Tue, Oct 26, 2021 at 10:17 AM Pierre
                    Seize <<a href="mailto:pierre.seize@onera.fr"
                      target="_blank" moz-do-not-send="true">pierre.seize@onera.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 bgcolor="#FFFFFF">
                        <p>Hi, I had the idea to try and renumber my
                          mesh cells, as I've heard it's better:
                          "neighbouring cells are stored next to one
                          another, and memory access are faster".</p>
                        Right now, I load the mesh then I distribute it
                        over the processes. I thought I'd try to permute
                        the numbering between those two steps :<br>
                        <blockquote><tt>DMPlexCreateFromFile</tt><tt><br>
                          </tt><tt>DMPlexGetOrdering</tt><tt><br>
                          </tt><tt>DMPlexPermute</tt><tt><br>
                          </tt><tt>DMPlexDistribute</tt></blockquote>
                        <p>but that gives me an error when it runs on
                          more than one process:</p>
                        <p><tt>[0]PETSC ERROR: ---------------------
                            Error Message
--------------------------------------------------------------</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: No support for this
                            operation for this object type</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: Number of dofs for
                            point 0 in the local section should be
                            positive</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: See <a
                              href="https://petsc.org/release/faq/"
                              target="_blank" moz-do-not-send="true">https://petsc.org/release/faq/</a>
                            for trouble shooting.</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: Petsc Release Version
                            3.16.0, unknown </tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: ./build/bin/yanss on
                            a  named ldmpe202z.onera by pseize Tue Oct
                            26 16:03:33 2021</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: Configure options
                            --PETSC_ARCH=arch-ld-gcc --download-metis
                            --download-parmetis --prefix=~/.local
                            --with-cgns</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: #1
                            PetscPartitionerDMPlexPartition() at
                            /stck/pseize/softwares/petsc/src/dm/impls/plex/plexpartition.c:720</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: #2 DMPlexDistribute()
                            at
                            /stck/pseize/softwares/petsc/src/dm/impls/plex/plexdistribute.c:1630</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: #3 MeshLoadFromFile()
                            at src/spatial.c:689</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: #4 main() at
                            src/main.c:22</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: PETSc Option Table
                            entries:</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: -draw_comp 0</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: -mesh data/box.msh</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: -mesh_view draw</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: -riemann anrs</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: -ts_max_steps 100</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: -vec_view_partition</tt><tt><br>
                          </tt><tt>[0]PETSC ERROR: ----------------End
                            of Error Message -------send entire error
                            message to <a
                              href="mailto:petsc-maint@mcs.anl.gov"
                              target="_blank" moz-do-not-send="true">petsc-maint@mcs.anl.gov</a>----------</tt><br>
                          <br>
                        </p>
                        <p>I checked and before I tried to reorder the
                          mesh, the <tt>dm-></tt><tt>localSection</tt>
                          was <tt>NULL</tt> before entering <tt>DMPlexDistribute</tt>,
                          and I was able to fix the error with <tt>DMSetLocalSection(dm,
                            NULL)</tt> after <tt>DMPlexPermute</tt>,
                          but it doesn't seems it's the right way to do
                          what I want. Does someone have any advice ?</p>
                      </div>
                    </blockquote>
                    <div>Oh, this is probably me trying to be too
                      clever. If a local section is defined, then I try
                      to use the number of dofs in it to load balance
                      better.</div>
                    <div>There should never be a negative number of dofs
                      in the local section (a global section uses this
                      to indicate a dof owned by another process).</div>
                    <div>So eliminating the local section will
                      definitely fix that error.</div>
                    <div><br>
                    </div>
                    <div>Now the question of how you got a local
                      section. DMPlexPermute()  does not create one, so
                      it seems like you had one ahead of time, and that</div>
                    <div>the values were not valid.</div>
                  </div>
                </div>
              </blockquote>
              <br>
              <tt>DMPlexPermute</tt> calls <tt>DMGetLocalSection</tt>,
              which creates <tt>dm->localSection</tt> if it's <tt>NULL</tt>,
              so before <tt>DMPlexPermute</tt> my <tt>dm->localSection</tt>
              is <tt>NULL</tt>, and after it is set. Because of that I
              enter the if in <tt>src/dm/impls/plex/plexpartition.c:707</tt>
              and then I got the error.<br>
              <tt> </tt>If i have a "wrong" <tt>dm->localSection</tt>,
              I think it has to come from <tt>DMPlexPermute</tt>.<br>
              <br>
              <blockquote type="cite">
                <div dir="ltr">
                  <div class="gmail_quote">
                    <div>Note that you can probably get rid of some of
                      the loading code using</div>
                    <div><br>
                    </div>
                    <div>  DMCreate(comm, &dm);</div>
                    <div>  DMSetType(dm, DMPLEX);</div>
                    <div>  DMSetFromOptions(dm);</div>
                    <div>  DMViewFromOptions(dm, NULL, "-mesh_view");</div>
                    <div><br>
                    </div>
                    <div>and use</div>
                    <div><br>
                    </div>
                    <div>  -dm_plex_filename databox,msh -mesh_view</div>
                  </div>
                </div>
              </blockquote>
              <br>
              My loading code is already small, but just to make sure I
              wrote this minimal example:<br>
              <tt><br>
                int main(int argc, char **argv){<br>
                  PetscErrorCode ierr;<br>
                <br>
                  ierr = PetscInitialize(&argc, &argv, NULL,
                help); if (ierr) return ierr;<br>
                <br>
                  DM dm, foo_dm;<br>
                  ierr = DMCreate(PETSC_COMM_WORLD, &dm);
                CHKERRQ(ierr);<br>
                  ierr = DMSetType(dm, DMPLEX); CHKERRQ(ierr);<br>
                  ierr = DMSetFromOptions(dm); CHKERRQ(ierr);<br>
                <br>
                  IS perm;<br>
                  ierr = DMPlexGetOrdering(dm, NULL, NULL, &perm);
                CHKERRQ(ierr);<br>
                  ierr = DMPlexPermute(dm, perm, &foo_dm);
                CHKERRQ(ierr);<br>
                  if (foo_dm) {<br>
                    ierr = DMDestroy(&dm); CHKERRQ(ierr);<br>
                    dm = foo_dm;<br>
                  }<br>
                  ierr = DMPlexDistribute(dm, 2, NULL, &foo_dm);
                CHKERRQ(ierr);<br>
                  if (foo_dm) {<br>
                    ierr = DMDestroy(&dm); CHKERRQ(ierr);<br>
                    dm = foo_dm;<br>
                  }<br>
                <br>
                  ierr = ISDestroy(&perm); CHKERRQ(ierr);<br>
                  ierr = DMDestroy(&dm); CHKERRQ(ierr);<br>
                  ierr = PetscFinalize();<br>
                  return ierr;<br>
                }</tt><br>
              <br>
              ran with <tt>mpiexec -n 2 ./build/bin/yanss
                -dm_plex_filename data/box.msh</tt>. The mesh is a 2D
              box from GMSH but I've got the same result with any mesh
              I've tried. It runs fine with 1 process but gives the
              previous error for more processes.<br>
            </div>
          </blockquote>
          <div><br>
          </div>
          <div>Hi Pierre,</div>
          <div><br>
          </div>
          <div>You are right. This is my bug. Here is the fix:</div>
          <div><br>
          </div>
          <div>  <a
              href="https://gitlab.com/petsc/petsc/-/merge_requests/4504"
              moz-do-not-send="true">https://gitlab.com/petsc/petsc/-/merge_requests/4504</a></div>
          <div><br>
          </div>
          <div>Is it possible to try this branch?</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 bgcolor="#FFFFFF"> Pierre<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/%7Eknepley/"
                        target="_blank" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                    </div>
                  </div>
                </div>
              </div>
            </div>
          </div>
        </div>
      </div>
    </blockquote>
    <br>
  </body>
</html>