<div dir="ltr"><div dir="ltr">On Wed, Aug 7, 2019 at 11:32 AM Lawrence Mitchell <<a href="mailto:wence@gmx.li">wence@gmx.li</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">On Wed, 7 Aug 2019 at 13:52, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<br>
><br>
> On Wed, Aug 7, 2019 at 7:13 AM Lawrence Mitchell via petsc-dev <<a href="mailto:petsc-dev@mcs.anl.gov" target="_blank">petsc-dev@mcs.anl.gov</a>> wrote:<br>
>><br>
>> Dear petsc-dev,<br>
>><br>
>> I would like to run with a geometric multigrid hierarchy built out of DMPlex objects. On the coarse grids, I would like to reduce the size of the communicator being used. Since I build the hierarchy by regularly refining the coarse grid problem, this means I need a way of redistributing a DMPlex object from commA to commB.<br>
>><br>
>> (I note that DMRefine has signature DMRefine(oldDM, newcomm, newDM), but as far as I can tell no concrete implementation of DMRefine has ever taken any notice of newcomm, and the newDM always ends up on the same communicator as the oldDM).<br>
>><br>
>> OK, so I think what I want to use is DMPlexDistribute. This does many-to-many redistribution of DMPlex objects, but has a rather constrained interface:<br>
>><br>
>> DMPlexDistribute(oldDM, overlap, *migrationSF, *newDM)<br>
>><br>
>> It works by taking the oldDM, assuming that the communicator of that thing defines the number of partitions I want, partitioning the mesh dual graph, and migrating onto a newDM, with the same number of processes as the oldDM.<br>
>><br>
>> So what are my options if I want to distribute from oldDM onto newDM where the communicators of the oldDM and the putative target newDM /don't/ match?<br>
>><br>
>> Given my use case above, right now I can assume that MPI_Group_translate_ranks(oldcomm, newcomm) will work (although I can think of future use cases where it doesn't).<br>
>><br>
>> In an ideal world, I think the interface for DMPlexDistribute should be something like:<br>
>><br>
>> DMPlexDistribute(oldDM, overlap, newComm, *migrationSF, *newDM)<br>
>><br>
>> where oldDM->comm and newComm are possibly disjoint sets of ranks, and the call is collective over the union of the groups of the two communicators.<br>
>><br>
>> This would work by then asking the partitioner to construct a partition based on the size of newComm (rather than oldDM->comm) and doing all of the migration.<br>
>><br>
>> I started to look at this, and swiftly threw my hands in the air, because the implicit assumption of the communicators matching is everywhere.<br>
>><br>
>> So my stopgap solution would be (given that I know oldDM->comm \subset newComm) is to take my oldDM and move it onto newComm before calling the existing DMPlexDistribute.<br>
>><br>
>> Does this sound like a reasonable approach? Effectively I need to create "empty" DM objects on all of the participating ranks in newComm that are not in oldComm: what information do they need to contain to align with oldDM?<br>
>><br>
>><br>
>> Does anyone else just have code for doing this lying around? Am I completely barking up the wrong tree?<br>
><br>
><br>
> That is how you do it. We are solidifying this pattern, as you can see from Junchao's new example for pushing a Vec onto a subcomm.<br>
<br>
OK, so I think I am getting there. Presently I am abusing<br>
DMPlexCreateFromDAG to migrate a DM on oldComm onto newComm, but this<br>
is very fragile. I attach what I have right now. You have to run it<br>
with PTSCOTCH, because parmetis refuses to partition graphs with no<br>
vertices on a process.: again, this would be avoided if the<br>
partitioning was done on a source communicator with a number of<br>
partitions given by a target communicator, anyway.<br></blockquote><div><br></div><div>Sorry I am just getting to this. Its a great example. I was thinking of just pushing this stuff in the library, but</div><div>I had the following thought. What if we reused DMClone() to stick things on another Comm, since we do</div><div>not actually want a copy. The internals would then have to contend with a NULL impl, which might be a lot</div><div>of trouble. I was going to try it out in a branch. It seems more elegant to me.</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">
This example builds a simple DMPlex on a subcommunicator (built with<br>
processes from COMM_WORLD where rank % 3 == 0), and then refines it,<br>
and redistributes it onto COMM_WORLD.<br>
<br>
Running with:<br>
<br>
$ mpiexec -n 4 ./dmplex-redist -input_dm_view -refined_dm_view<br>
-copied_dm_view -redistributed_dm_view<br>
<br>
Builds a box mesh on 1 rank, redistributes it onto 2, refines that DM,<br>
then expands it to live on 4 ranks, and subsequently redistributes<br>
that.<br>
<br>
DM Object: 2 MPI processes<br>
  type: plex<br>
DM_0x84000007_0 in 3 dimensions:<br>
  0-cells: 27 0<br>
  1-cells: 54 0<br>
  2-cells: 36 0<br>
  3-cells: 8 0<br>
Labels:<br>
  depth: 4 strata with value/size (0 (27), 1 (54), 2 (36), 3 (8))<br>
  Face Sets: 6 strata with value/size (6 (4), 5 (4), 3 (4), 4 (4), 1 (4), 2 (4))<br>
  marker: 1 strata with value/size (1 (72))<br>
DM Object: 2 MPI processes<br>
  type: plex<br>
DM_0x84000007_1 in 3 dimensions:<br>
  0-cells: 75 75<br>
  1-cells: 170 170<br>
  2-cells: 128 128<br>
  3-cells: 32 32<br>
Labels:<br>
  marker: 1 strata with value/size (1 (192))<br>
  Face Sets: 5 strata with value/size (1 (36), 3 (18), 4 (18), 5 (18), 6 (18))<br>
  depth: 4 strata with value/size (0 (75), 1 (170), 2 (128), 3 (32))<br>
DM Object: 4 MPI processes<br>
  type: plex<br>
DM_0xc4000003_0 in 3 dimensions:<br>
  0-cells: 75 0 0 75<br>
  1-cells: 170 0 0 170<br>
  2-cells: 128 0 0 128<br>
  3-cells: 32 0 0 32<br>
Labels:<br>
  depth: 4 strata with value/size (0 (75), 1 (170), 2 (128), 3 (32))<br>
Field PetscContainer_0xc4000003_1:<br>
  adjacency FEM<br>
DM Object: Parallel Mesh 4 MPI processes<br>
  type: plex<br>
Parallel Mesh in 3 dimensions:<br>
  0-cells: 45 45 45 45<br>
  1-cells: 96 96 96 96<br>
  2-cells: 68 68 68 68<br>
  3-cells: 16 16 16 16<br>
Labels:<br>
  depth: 4 strata with value/size (0 (45), 1 (96), 2 (68), 3 (16))<br>
Field PetscContainer_0xc4000003_2:<br>
  adjacency FEM<br>
<br>
I would ideally like some help in figuring out a better way of<br>
"expanding" the DM onto the new communicator (rather than using<br>
DMPlexCreateFromDAG as I do now).<br>
<br>
I think one wants to think about the interface for these<br>
redistribution functions in general. It seems that they want to<br>
broadly match MPI_Intercomm_create. So something like:<br>
<br>
DMPlexDistribute(oldDM, peerCommunicator, newCommunicator, overlap,<br>
*migrationSF, *newDM)<br>
<br>
This is collective over peerCommunicator and returns a newDM on<br>
newCommunicator, and a migrationSF on peerCommunicator.<br>
<br>
Note that I do not think it is safe to assume that PETSC_COMM_WORLD is<br>
always suitable as peerCommunicator.<br>
<br>
> I think the right way to do this would be to implement the hooks in PCTELESCOPE for DMPlex. Dave and I have talked about this and<br>
> it should be exactly the same work as you propose above, but it would allow you to use the command line, do this recursively, interact nicely<br>
> with the solvers, etc. I can help.<br>
<br>
The telescope side of things now exists (c.f. 8d9f7141f511) to some<br>
degree. But to do that, one needs the redistribution.<br>
<br>
Cheers,<br>
<br>
Lawrence<br>
</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>