<div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr">Alright, I think the problem had to do with keeping track of what DM you were looking at. This code increases the overlap of an initial DM:<div><br></div><div>static char help[] = "Tests plex distribution and overlaps.\n";<br><br>#include <petsc/private/dmpleximpl.h><br><br>int main (int argc, char * argv[]) {<br><br>  DM             dm, dm2;<br>  PetscInt       overlap;<br>  MPI_Comm       comm;<br>  PetscErrorCode ierr;<br><br>  ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;<br>  comm = PETSC_COMM_WORLD;<br><br>  ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);<br>  ierr = DMSetFromOptions(dm);CHKERRQ(ierr);<br>  ierr = PetscObjectSetName((PetscObject) dm, "Initial DM");CHKERRQ(ierr);<br>  ierr = DMViewFromOptions(dm, NULL, "-initial_dm_view");CHKERRQ(ierr);<br><br>  ierr = DMPlexGetOverlap(dm, &overlap);CHKERRQ(ierr);<br>  ierr = DMPlexDistributeOverlap(dm, overlap+1, NULL, &dm2);CHKERRQ(ierr);<br>  ierr = PetscObjectSetName((PetscObject) dm2, "More Overlap DM");CHKERRQ(ierr);<br>  ierr = DMViewFromOptions(dm2, NULL, "-over_dm_view");CHKERRQ(ierr);<br><br>  ierr = DMDestroy(&dm2);CHKERRQ(ierr);<br>  ierr = DMDestroy(&dm);CHKERRQ(ierr);<br>  ierr = PetscFinalize();<br>  return ierr;<br>}<br></div><div><br></div><div>and when we run it we get the expected result</div><div><br></div><div>master *:~/Downloads/tmp/Nicolas$ /PETSc3/petsc/apple/bin/mpiexec -n 2 ./test_overlap -initial_dm_view -dm_plex_box_faces 5,5 -dm_distribute -dm_distribute_overlap 1 -over_dm_view<br>DM Object: Initial DM 2 MPI processes<br>  type: plex<br>Initial DM in 2 dimensions:<br>  0-cells: 29 29<br>  1-cells: 65 65<br>  2-cells: 37 37<br>Labels:<br>  depth: 3 strata with value/size (0 (29), 1 (65), 2 (37))<br>  celltype: 3 strata with value/size (0 (29), 1 (65), 3 (37))<br>  marker: 1 strata with value/size (1 (27))<br>  Face Sets: 1 strata with value/size (1 (13))<br>DM Object: More Overlap DM 2 MPI processes<br>  type: plex<br>More Overlap DM in 2 dimensions:<br>  0-cells: 36 36<br>  1-cells: 85 85<br>  2-cells: 50 50<br>Labels:<br>  depth: 3 strata with value/size (0 (36), 1 (85), 2 (50))<br>  celltype: 3 strata with value/size (0 (36), 1 (85), 3 (50))<br>  marker: 1 strata with value/size (1 (40))<br>  Face Sets: 1 strata with value/size (1 (20))<br></div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div></div></div></div></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Mar 31, 2021 at 12:57 PM Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@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">Okay, let me show a really simple example that gives the expected result before I figure out what is going wrong for you. This code<div><br></div><div>static char help[] = "Tests plex distribution and overlaps.\n";<br><br>#include <petsc/private/dmpleximpl.h><br><br>int main (int argc, char * argv[]) {<br></div><div>  DM                    dm;</div><div>  MPI_Comm       comm;</div><div>  PetscErrorCode ierr;</div><div><br></div><div>  ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;<br>  comm = PETSC_COMM_WORLD;<br><br>  ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);<br>  ierr = DMSetFromOptions(dm);CHKERRQ(ierr);<br>  ierr = PetscObjectSetName((PetscObject) dm, "Initial DM");CHKERRQ(ierr);<br>  ierr = DMViewFromOptions(dm, NULL, "-initial_dm_view");CHKERRQ(ierr);<br></div><div>  ierr = DMDestroy(&dm);CHKERRQ(ierr);<br>  ierr = PetscFinalize();<br>  return ierr;<br>}<br></div><div><br></div><div>can do all the overlap tests. For example, you can run it naively and get a serial mesh</div><div><br></div><div>master *:~/Downloads/tmp/Nicolas$ /PETSc3/petsc/apple/bin/mpiexec -n 2 ./test_overlap -initial_dm_view -dm_plex_box_faces 5,5<br>DM Object: Initial DM 2 MPI processes<br>  type: plex<br>Initial DM in 2 dimensions:<br>  0-cells: 36 0<br>  1-cells: 85 0<br>  2-cells: 50 0<br>Labels:<br>  celltype: 3 strata with value/size (0 (36), 3 (50), 1 (85))<br>  depth: 3 strata with value/size (0 (36), 1 (85), 2 (50))<br>  marker: 1 strata with value/size (1 (40))<br>  Face Sets: 1 strata with value/size (1 (20))</div><div><br></div><div>Then run it telling Plex to distribute after creating the mesh</div><div><br>master *:~/Downloads/tmp/Nicolas$ /PETSc3/petsc/apple/bin/mpiexec -n 2 ./test_overlap -initial_dm_view -dm_plex_box_faces 5,5 -dm_distribute<br>DM Object: Initial DM 2 MPI processes<br>  type: plex<br>Initial DM in 2 dimensions:<br>  0-cells: 21 21<br>  1-cells: 45 45<br>  2-cells: 25 25<br>Labels:<br>  depth: 3 strata with value/size (0 (21), 1 (45), 2 (25))<br>  celltype: 3 strata with value/size (0 (21), 1 (45), 3 (25))<br>  marker: 1 strata with value/size (1 (21))<br>  Face Sets: 1 strata with value/size (1 (10))</div><div><br></div><div>The get the same thing back with overlap = 0</div><div><br>master *:~/Downloads/tmp/Nicolas$ /PETSc3/petsc/apple/bin/mpiexec -n 2 ./test_overlap -initial_dm_view -dm_plex_box_faces 5,5 -dm_distribute -dm_distribute_overlap 0<br>DM Object: Initial DM 2 MPI processes<br>  type: plex<br>Initial DM in 2 dimensions:<br>  0-cells: 21 21<br>  1-cells: 45 45<br>  2-cells: 25 25<br>Labels:<br>  depth: 3 strata with value/size (0 (21), 1 (45), 2 (25))<br>  celltype: 3 strata with value/size (0 (21), 1 (45), 3 (25))<br>  marker: 1 strata with value/size (1 (21))<br>  Face Sets: 1 strata with value/size (1 (10))</div><div><br></div><div>and get larger local meshes with overlap = 1</div><div><br>master *:~/Downloads/tmp/Nicolas$ /PETSc3/petsc/apple/bin/mpiexec -n 2 ./test_overlap -initial_dm_view -dm_plex_box_faces 5,5 -dm_distribute -dm_distribute_overlap 1<br>DM Object: Initial DM 2 MPI processes<br>  type: plex<br>Initial DM in 2 dimensions:<br>  0-cells: 29 29<br>  1-cells: 65 65<br>  2-cells: 37 37<br>Labels:<br>  depth: 3 strata with value/size (0 (29), 1 (65), 2 (37))<br>  celltype: 3 strata with value/size (0 (29), 1 (65), 3 (37))<br>  marker: 1 strata with value/size (1 (27))<br>  Face Sets: 1 strata with value/size (1 (13))<br></div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Mar 31, 2021 at 12:22 PM Nicolas Barral <<a href="mailto:nicolas.barral@math.u-bordeaux.fr" target="_blank">nicolas.barral@math.u-bordeaux.fr</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>
@+<br>
<br>
-- <br>
Nicolas<br>
<br>
On 31/03/2021 17:51, Matthew Knepley wrote:<br>
> On Sat, Mar 27, 2021 at 9:27 AM Nicolas Barral <br>
> <<a href="mailto:nicolas.barral@math.u-bordeaux.fr" target="_blank">nicolas.barral@math.u-bordeaux.fr</a> <br>
> <mailto:<a href="mailto:nicolas.barral@math.u-bordeaux.fr" target="_blank">nicolas.barral@math.u-bordeaux.fr</a>>> wrote:<br>
> <br>
>     Hi all,<br>
> <br>
>     First, I'm not sure I understand what the overlap parameter in<br>
>     DMPlexDistributeOverlap does. I tried the following: generate a small<br>
>     mesh on 1 rank with DMPlexCreateBoxMesh, then distribute it with<br>
>     DMPlexDistribute. At this point I have two nice partitions, with shared<br>
>     vertices and no overlapping cells. Then I call DMPlexDistributeOverlap<br>
>     with the overlap parameter set to 0 or 1, and get the same resulting<br>
>     plex in both cases. Why is that ?<br>
> <br>
> <br>
> The overlap parameter says how many cell adjacencies to go out. You <br>
> should not get the same<br>
> mesh out. We have lots of examples that use this. If you send your small <br>
> example, I can probably<br>
> tell you what is happening.<br>
> <br>
<br>
Ok so I do have a small example on that and the DMClone thing I set up <br>
to understand! I attach it to the email.<br>
<br>
For the overlap, you can change the overlap constant at the top of the <br>
file. With OVERLAP=0 or 1, the distributed overlapping mesh (shown using <br>
-over_dm_view, it's DMover) are the same, and different from the mesh <br>
before distributing the overlap (shown using -distrib_dm_view). For <br>
larger overlap values they're different.<br>
<br>
The process is:<br>
1/ create a DM dm on 1 rank<br>
2/ clone dm into dm2<br>
3/ distribute dm<br>
4/ clone dm into dm3<br>
5/ distribute dm overlap<br>
<br>
I print all the DMs after each step. dm has a distributed overlap, dm2 <br>
is not distributed, dm3 is distributed but without overlap. Since <br>
distribute and distributeOverlap create new DMs, I don't seem have a <br>
problem with the shallow copies.<br>
<br>
<br>
>     Second, I'm wondering what would be a good way to handle two overlaps<br>
>     and associated local vectors. In my adaptation code, the remeshing<br>
>     library requires a non-overlapping mesh, while the refinement criterion<br>
>     computation is based on hessian computations, which require a layer of<br>
>     overlap. What I can do is clone the dm before distributing the overlap,<br>
>     then manage two independent plex objects with their own local sections<br>
>     etc. and copy/trim local vectors manually. Is there a more automatic<br>
>     way<br>
>     to do this ?<br>
> <br>
> <br>
> DMClone() is a shallow copy, so that will not work. You would maintain <br>
> two different Plexes, overlapping<br>
> and non-overlapping, with their own sections and vecs. Are you sure you <br>
> need to keep around the non-overlapping one?<br>
> Maybe if I understood what operations you want to work, I could say <br>
> something more definitive.<br>
> <br>
I need to be able to pass the non-overlapping mesh to the remesher. I <br>
can either maintain 2 plexes, or trim the overlapping plex when I create <br>
the arrays I give to the remesher. I'm not sure which is the best/worst ?<br>
<br>
Thanks<br>
<br>
-- <br>
Nicolas<br>
<br>
<br>
>    Thanks,<br>
> <br>
>       Matt<br>
> <br>
>     Thanks<br>
> <br>
>     -- <br>
>     Nicolas<br>
> <br>
> <br>
> <br>
> -- <br>
> What most experimenters take for granted before they begin their <br>
> experiments is infinitely more interesting than any results to which <br>
> their experiments lead.<br>
> -- Norbert Wiener<br>
> <br>
> <a href="https://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~knepley/</a> <<a href="http://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">http://www.cse.buffalo.edu/~knepley/</a>><br>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><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>
</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>