[petsc-dev] DMGetNeighbors and DMPlex

Pierre Jolivet pierre.jolivet at enseeiht.fr
Mon Apr 13 05:13:48 CDT 2020


Hello,
I’m trying to figure out why the result of DMGetNeighbors when a DMPlex is distributed with no overlap can be nonsymmetric.
Is this intended? If so, is there some other easy way to get the equivalent list as in the overlapping case (I’m interested in the list of ranks sharing at least one point)?
Here is a MWE.
$ patch -p1 < patch.txt
$ cd src/dm/impls/plex/tests/
$ make ex1
$ mpirun -n 2 ./ex1 -dim 2 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3 -test_redistribute
1 <— how come rank 0 has 1 neighbor
0 <— but rank 1 has 0 neighbor?
$ mpirun -n 2 ./ex1 -dim 2 -cell_simplex 0 -interpolate -domain_shape box -domain_box_sizes 3,3 -overlap 1
1
1
Overlap: 1

Thanks,
Pierre

diff --git a/src/dm/impls/plex/tests/ex1.c b/src/dm/impls/plex/tests/ex1.c
index 4a5d051a66..ac2ce8eb20 100644
--- a/src/dm/impls/plex/tests/ex1.c
+++ b/src/dm/impls/plex/tests/ex1.c
@@ -3,0 +4 @@ static char help[] = "Tests various DMPlex routines to construct, refine and dis
+#include <petscsf.h>
@@ -420,0 +422,8 @@ PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
+      PetscInt nranks;
+      const PetscInt* ranks;
+      PetscSF sf;
+      ierr = DMGetPointSF(distributedMesh, &sf);CHKERRQ(ierr);
+      ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
+      ierr = DMGetNeighbors(distributedMesh, &nranks, &ranks);CHKERRQ(ierr);
+      ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%D\n",nranks);CHKERRQ(ierr);
+      ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
@@ -434,0 +444,8 @@ PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
+      PetscInt nranks;
+      const PetscInt* ranks;
+      PetscSF sf;
+      ierr = DMGetPointSF(overlapMesh, &sf);CHKERRQ(ierr);
+      ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
+      ierr = DMGetNeighbors(overlapMesh, &nranks, &ranks);CHKERRQ(ierr);
+      ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%D\n",nranks);CHKERRQ(ierr);
+      ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);


More information about the petsc-dev mailing list