1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 370034214SMatthew G. Knepley 43c1f0c11SLawrence Mitchell /*@C 53c1f0c11SLawrence Mitchell DMPlexSetAdjacencyUser - Define adjacency in the mesh using a user-provided callback 63c1f0c11SLawrence Mitchell 73c1f0c11SLawrence Mitchell Input Parameters: 83c1f0c11SLawrence Mitchell + dm - The DM object 93c1f0c11SLawrence Mitchell . user - The user callback, may be NULL (to clear the callback) 103c1f0c11SLawrence Mitchell - ctx - context for callback evaluation, may be NULL 113c1f0c11SLawrence Mitchell 123c1f0c11SLawrence Mitchell Level: advanced 133c1f0c11SLawrence Mitchell 143c1f0c11SLawrence Mitchell Notes: 153c1f0c11SLawrence Mitchell The caller of DMPlexGetAdjacency may need to arrange that a large enough array is available for the adjacency. 163c1f0c11SLawrence Mitchell 173c1f0c11SLawrence Mitchell Any setting here overrides other configuration of DMPlex adjacency determination. 183c1f0c11SLawrence Mitchell 193c1f0c11SLawrence Mitchell .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexGetAdjacency(), DMPlexGetAdjacencyUser() 203c1f0c11SLawrence Mitchell @*/ 213c1f0c11SLawrence Mitchell PetscErrorCode DMPlexSetAdjacencyUser(DM dm,PetscErrorCode (*user)(DM,PetscInt,PetscInt*,PetscInt[],void*),void *ctx) 223c1f0c11SLawrence Mitchell { 233c1f0c11SLawrence Mitchell DM_Plex *mesh = (DM_Plex *)dm->data; 243c1f0c11SLawrence Mitchell 253c1f0c11SLawrence Mitchell PetscFunctionBegin; 263c1f0c11SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 273c1f0c11SLawrence Mitchell mesh->useradjacency = user; 283c1f0c11SLawrence Mitchell mesh->useradjacencyctx = ctx; 293c1f0c11SLawrence Mitchell PetscFunctionReturn(0); 303c1f0c11SLawrence Mitchell } 313c1f0c11SLawrence Mitchell 323c1f0c11SLawrence Mitchell /*@C 333c1f0c11SLawrence Mitchell DMPlexGetAdjacencyUser - get the user-defined adjacency callback 343c1f0c11SLawrence Mitchell 353c1f0c11SLawrence Mitchell Input Parameter: 363c1f0c11SLawrence Mitchell . dm - The DM object 373c1f0c11SLawrence Mitchell 383c1f0c11SLawrence Mitchell Output Parameters: 393c1f0c11SLawrence Mitchell - user - The user callback 403c1f0c11SLawrence Mitchell - ctx - context for callback evaluation 413c1f0c11SLawrence Mitchell 423c1f0c11SLawrence Mitchell Level: advanced 433c1f0c11SLawrence Mitchell 443c1f0c11SLawrence Mitchell .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexGetAdjacency(), DMPlexSetAdjacencyUser() 453c1f0c11SLawrence Mitchell @*/ 463c1f0c11SLawrence Mitchell PetscErrorCode DMPlexGetAdjacencyUser(DM dm, PetscErrorCode (**user)(DM,PetscInt,PetscInt*,PetscInt[],void*), void **ctx) 473c1f0c11SLawrence Mitchell { 483c1f0c11SLawrence Mitchell DM_Plex *mesh = (DM_Plex *)dm->data; 493c1f0c11SLawrence Mitchell 503c1f0c11SLawrence Mitchell PetscFunctionBegin; 513c1f0c11SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 523c1f0c11SLawrence Mitchell if (user) *user = mesh->useradjacency; 533c1f0c11SLawrence Mitchell if (ctx) *ctx = mesh->useradjacencyctx; 543c1f0c11SLawrence Mitchell PetscFunctionReturn(0); 553c1f0c11SLawrence Mitchell } 563c1f0c11SLawrence Mitchell 5770034214SMatthew G. Knepley /*@ 5870034214SMatthew G. Knepley DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first 5970034214SMatthew G. Knepley 6070034214SMatthew G. Knepley Input Parameters: 6170034214SMatthew G. Knepley + dm - The DM object 6270034214SMatthew G. Knepley - useCone - Flag to use the cone first 6370034214SMatthew G. Knepley 6470034214SMatthew G. Knepley Level: intermediate 6570034214SMatthew G. Knepley 6670034214SMatthew G. Knepley Notes: 6770034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 684b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 6970034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 7070034214SMatthew G. Knepley 7170034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 7270034214SMatthew G. Knepley @*/ 7370034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone) 7470034214SMatthew G. Knepley { 7570034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 7670034214SMatthew G. Knepley 7770034214SMatthew G. Knepley PetscFunctionBegin; 7870034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7970034214SMatthew G. Knepley mesh->useCone = useCone; 8070034214SMatthew G. Knepley PetscFunctionReturn(0); 8170034214SMatthew G. Knepley } 8270034214SMatthew G. Knepley 8370034214SMatthew G. Knepley /*@ 8470034214SMatthew G. Knepley DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first 8570034214SMatthew G. Knepley 8670034214SMatthew G. Knepley Input Parameter: 8770034214SMatthew G. Knepley . dm - The DM object 8870034214SMatthew G. Knepley 8970034214SMatthew G. Knepley Output Parameter: 9070034214SMatthew G. Knepley . useCone - Flag to use the cone first 9170034214SMatthew G. Knepley 9270034214SMatthew G. Knepley Level: intermediate 9370034214SMatthew G. Knepley 9470034214SMatthew G. Knepley Notes: 9570034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 964b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 9770034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 9870034214SMatthew G. Knepley 9970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 10070034214SMatthew G. Knepley @*/ 10170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone) 10270034214SMatthew G. Knepley { 10370034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 10470034214SMatthew G. Knepley 10570034214SMatthew G. Knepley PetscFunctionBegin; 10670034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10770034214SMatthew G. Knepley PetscValidIntPointer(useCone, 2); 10870034214SMatthew G. Knepley *useCone = mesh->useCone; 10970034214SMatthew G. Knepley PetscFunctionReturn(0); 11070034214SMatthew G. Knepley } 11170034214SMatthew G. Knepley 11270034214SMatthew G. Knepley /*@ 11370034214SMatthew G. Knepley DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure 11470034214SMatthew G. Knepley 11570034214SMatthew G. Knepley Input Parameters: 11670034214SMatthew G. Knepley + dm - The DM object 11770034214SMatthew G. Knepley - useClosure - Flag to use the closure 11870034214SMatthew G. Knepley 11970034214SMatthew G. Knepley Level: intermediate 12070034214SMatthew G. Knepley 12170034214SMatthew G. Knepley Notes: 12270034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 1234b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 12470034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 12570034214SMatthew G. Knepley 12670034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 12770034214SMatthew G. Knepley @*/ 12870034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure) 12970034214SMatthew G. Knepley { 13070034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 13170034214SMatthew G. Knepley 13270034214SMatthew G. Knepley PetscFunctionBegin; 13370034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13470034214SMatthew G. Knepley mesh->useClosure = useClosure; 13570034214SMatthew G. Knepley PetscFunctionReturn(0); 13670034214SMatthew G. Knepley } 13770034214SMatthew G. Knepley 13870034214SMatthew G. Knepley /*@ 13970034214SMatthew G. Knepley DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure 14070034214SMatthew G. Knepley 14170034214SMatthew G. Knepley Input Parameter: 14270034214SMatthew G. Knepley . dm - The DM object 14370034214SMatthew G. Knepley 14470034214SMatthew G. Knepley Output Parameter: 14570034214SMatthew G. Knepley . useClosure - Flag to use the closure 14670034214SMatthew G. Knepley 14770034214SMatthew G. Knepley Level: intermediate 14870034214SMatthew G. Knepley 14970034214SMatthew G. Knepley Notes: 15070034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 1514b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 15270034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 15370034214SMatthew G. Knepley 15470034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 15570034214SMatthew G. Knepley @*/ 15670034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure) 15770034214SMatthew G. Knepley { 15870034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 15970034214SMatthew G. Knepley 16070034214SMatthew G. Knepley PetscFunctionBegin; 16170034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 16270034214SMatthew G. Knepley PetscValidIntPointer(useClosure, 2); 16370034214SMatthew G. Knepley *useClosure = mesh->useClosure; 16470034214SMatthew G. Knepley PetscFunctionReturn(0); 16570034214SMatthew G. Knepley } 16670034214SMatthew G. Knepley 1678b0b4c70SToby Isaac /*@ 168a17985deSToby Isaac DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints. 1698b0b4c70SToby Isaac 1708b0b4c70SToby Isaac Input Parameters: 1718b0b4c70SToby Isaac + dm - The DM object 1725b317d89SToby Isaac - useAnchors - Flag to use the constraints. If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place. 1738b0b4c70SToby Isaac 1748b0b4c70SToby Isaac Level: intermediate 1758b0b4c70SToby Isaac 176a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors() 1778b0b4c70SToby Isaac @*/ 1785b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors) 1798b0b4c70SToby Isaac { 1808b0b4c70SToby Isaac DM_Plex *mesh = (DM_Plex *) dm->data; 1818b0b4c70SToby Isaac 1828b0b4c70SToby Isaac PetscFunctionBegin; 1838b0b4c70SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1845b317d89SToby Isaac mesh->useAnchors = useAnchors; 1858b0b4c70SToby Isaac PetscFunctionReturn(0); 1868b0b4c70SToby Isaac } 1878b0b4c70SToby Isaac 1888b0b4c70SToby Isaac /*@ 189a17985deSToby Isaac DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints. 1908b0b4c70SToby Isaac 1918b0b4c70SToby Isaac Input Parameter: 1928b0b4c70SToby Isaac . dm - The DM object 1938b0b4c70SToby Isaac 1948b0b4c70SToby Isaac Output Parameter: 1955b317d89SToby Isaac . useAnchors - Flag to use the closure. If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place. 1968b0b4c70SToby Isaac 1978b0b4c70SToby Isaac Level: intermediate 1988b0b4c70SToby Isaac 199a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors() 2008b0b4c70SToby Isaac @*/ 2015b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors) 2028b0b4c70SToby Isaac { 2038b0b4c70SToby Isaac DM_Plex *mesh = (DM_Plex *) dm->data; 2048b0b4c70SToby Isaac 2058b0b4c70SToby Isaac PetscFunctionBegin; 2068b0b4c70SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2075b317d89SToby Isaac PetscValidIntPointer(useAnchors, 2); 2085b317d89SToby Isaac *useAnchors = mesh->useAnchors; 2098b0b4c70SToby Isaac PetscFunctionReturn(0); 2108b0b4c70SToby Isaac } 2118b0b4c70SToby Isaac 21270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 21370034214SMatthew G. Knepley { 21470034214SMatthew G. Knepley const PetscInt *cone = NULL; 21570034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, coneSize, c; 21670034214SMatthew G. Knepley PetscErrorCode ierr; 21770034214SMatthew G. Knepley 21870034214SMatthew G. Knepley PetscFunctionBeginHot; 21970034214SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 22070034214SMatthew G. Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 2214b6b44bdSMatthew G. Knepley for (c = 0; c <= coneSize; ++c) { 2224b6b44bdSMatthew G. Knepley const PetscInt point = !c ? p : cone[c-1]; 22370034214SMatthew G. Knepley const PetscInt *support = NULL; 22470034214SMatthew G. Knepley PetscInt supportSize, s, q; 22570034214SMatthew G. Knepley 2264b6b44bdSMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 2274b6b44bdSMatthew G. Knepley ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 22870034214SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 229527e7439SSatish Balay for (q = 0; q < numAdj || ((void)(adj[numAdj++] = support[s]),0); ++q) { 23070034214SMatthew G. Knepley if (support[s] == adj[q]) break; 23170034214SMatthew G. Knepley } 23270034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 23370034214SMatthew G. Knepley } 23470034214SMatthew G. Knepley } 23570034214SMatthew G. Knepley *adjSize = numAdj; 23670034214SMatthew G. Knepley PetscFunctionReturn(0); 23770034214SMatthew G. Knepley } 23870034214SMatthew G. Knepley 23970034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 24070034214SMatthew G. Knepley { 24170034214SMatthew G. Knepley const PetscInt *support = NULL; 24270034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, supportSize, s; 24370034214SMatthew G. Knepley PetscErrorCode ierr; 24470034214SMatthew G. Knepley 24570034214SMatthew G. Knepley PetscFunctionBeginHot; 24670034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 24770034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 2484b6b44bdSMatthew G. Knepley for (s = 0; s <= supportSize; ++s) { 2494b6b44bdSMatthew G. Knepley const PetscInt point = !s ? p : support[s-1]; 25070034214SMatthew G. Knepley const PetscInt *cone = NULL; 25170034214SMatthew G. Knepley PetscInt coneSize, c, q; 25270034214SMatthew G. Knepley 2534b6b44bdSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2544b6b44bdSMatthew G. Knepley ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 25570034214SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 256527e7439SSatish Balay for (q = 0; q < numAdj || ((void)(adj[numAdj++] = cone[c]),0); ++q) { 25770034214SMatthew G. Knepley if (cone[c] == adj[q]) break; 25870034214SMatthew G. Knepley } 25970034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 26070034214SMatthew G. Knepley } 26170034214SMatthew G. Knepley } 26270034214SMatthew G. Knepley *adjSize = numAdj; 26370034214SMatthew G. Knepley PetscFunctionReturn(0); 26470034214SMatthew G. Knepley } 26570034214SMatthew G. Knepley 26670034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[]) 26770034214SMatthew G. Knepley { 26870034214SMatthew G. Knepley PetscInt *star = NULL; 26970034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, starSize, s; 27070034214SMatthew G. Knepley PetscErrorCode ierr; 27170034214SMatthew G. Knepley 27270034214SMatthew G. Knepley PetscFunctionBeginHot; 27370034214SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 27470034214SMatthew G. Knepley for (s = 0; s < starSize*2; s += 2) { 27570034214SMatthew G. Knepley const PetscInt *closure = NULL; 27670034214SMatthew G. Knepley PetscInt closureSize, c, q; 27770034214SMatthew G. Knepley 27870034214SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 27970034214SMatthew G. Knepley for (c = 0; c < closureSize*2; c += 2) { 280527e7439SSatish Balay for (q = 0; q < numAdj || ((void)(adj[numAdj++] = closure[c]),0); ++q) { 28170034214SMatthew G. Knepley if (closure[c] == adj[q]) break; 28270034214SMatthew G. Knepley } 28370034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 28470034214SMatthew G. Knepley } 28570034214SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 28670034214SMatthew G. Knepley } 28770034214SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 28870034214SMatthew G. Knepley *adjSize = numAdj; 28970034214SMatthew G. Knepley PetscFunctionReturn(0); 29070034214SMatthew G. Knepley } 29170034214SMatthew G. Knepley 2925b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[]) 29370034214SMatthew G. Knepley { 29479bad331SMatthew G. Knepley static PetscInt asiz = 0; 2958b0b4c70SToby Isaac PetscInt maxAnchors = 1; 2968b0b4c70SToby Isaac PetscInt aStart = -1, aEnd = -1; 2978b0b4c70SToby Isaac PetscInt maxAdjSize; 2988b0b4c70SToby Isaac PetscSection aSec = NULL; 2998b0b4c70SToby Isaac IS aIS = NULL; 3008b0b4c70SToby Isaac const PetscInt *anchors; 3013c1f0c11SLawrence Mitchell DM_Plex *mesh = (DM_Plex *)dm->data; 30270034214SMatthew G. Knepley PetscErrorCode ierr; 30370034214SMatthew G. Knepley 30470034214SMatthew G. Knepley PetscFunctionBeginHot; 3055b317d89SToby Isaac if (useAnchors) { 306a17985deSToby Isaac ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr); 3078b0b4c70SToby Isaac if (aSec) { 3088b0b4c70SToby Isaac ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr); 30924c766afSToby Isaac maxAnchors = PetscMax(1,maxAnchors); 3108b0b4c70SToby Isaac ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr); 3118b0b4c70SToby Isaac ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 3128b0b4c70SToby Isaac } 3138b0b4c70SToby Isaac } 31479bad331SMatthew G. Knepley if (!*adj) { 31524c766afSToby Isaac PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd; 31679bad331SMatthew G. Knepley 31724c766afSToby Isaac ierr = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr); 31879bad331SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 31924c766afSToby Isaac ierr = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr); 32024c766afSToby Isaac coneSeries = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1; 32124c766afSToby Isaac supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1; 32224c766afSToby Isaac asiz = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries); 3238b0b4c70SToby Isaac asiz *= maxAnchors; 32424c766afSToby Isaac asiz = PetscMin(asiz,pEnd-pStart); 32579bad331SMatthew G. Knepley ierr = PetscMalloc1(asiz,adj);CHKERRQ(ierr); 32679bad331SMatthew G. Knepley } 32779bad331SMatthew G. Knepley if (*adjSize < 0) *adjSize = asiz; 3288b0b4c70SToby Isaac maxAdjSize = *adjSize; 3293c1f0c11SLawrence Mitchell if (mesh->useradjacency) { 3303c1f0c11SLawrence Mitchell ierr = mesh->useradjacency(dm, p, adjSize, *adj, mesh->useradjacencyctx);CHKERRQ(ierr); 3313c1f0c11SLawrence Mitchell } else if (useTransitiveClosure) { 33279bad331SMatthew G. Knepley ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr); 33370034214SMatthew G. Knepley } else if (useCone) { 33479bad331SMatthew G. Knepley ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr); 33570034214SMatthew G. Knepley } else { 33679bad331SMatthew G. Knepley ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr); 33770034214SMatthew G. Knepley } 3385b317d89SToby Isaac if (useAnchors && aSec) { 3398b0b4c70SToby Isaac PetscInt origSize = *adjSize; 3408b0b4c70SToby Isaac PetscInt numAdj = origSize; 3418b0b4c70SToby Isaac PetscInt i = 0, j; 3428b0b4c70SToby Isaac PetscInt *orig = *adj; 3438b0b4c70SToby Isaac 3448b0b4c70SToby Isaac while (i < origSize) { 3458b0b4c70SToby Isaac PetscInt p = orig[i]; 3468b0b4c70SToby Isaac PetscInt aDof = 0; 3478b0b4c70SToby Isaac 3488b0b4c70SToby Isaac if (p >= aStart && p < aEnd) { 3498b0b4c70SToby Isaac ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 3508b0b4c70SToby Isaac } 3518b0b4c70SToby Isaac if (aDof) { 3528b0b4c70SToby Isaac PetscInt aOff; 3538b0b4c70SToby Isaac PetscInt s, q; 3548b0b4c70SToby Isaac 3558b0b4c70SToby Isaac for (j = i + 1; j < numAdj; j++) { 3568b0b4c70SToby Isaac orig[j - 1] = orig[j]; 3578b0b4c70SToby Isaac } 3588b0b4c70SToby Isaac origSize--; 3598b0b4c70SToby Isaac numAdj--; 3608b0b4c70SToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 3618b0b4c70SToby Isaac for (s = 0; s < aDof; ++s) { 362527e7439SSatish Balay for (q = 0; q < numAdj || ((void)(orig[numAdj++] = anchors[aOff+s]),0); ++q) { 3638b0b4c70SToby Isaac if (anchors[aOff+s] == orig[q]) break; 3648b0b4c70SToby Isaac } 3658b0b4c70SToby Isaac if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 3668b0b4c70SToby Isaac } 3678b0b4c70SToby Isaac } 3688b0b4c70SToby Isaac else { 3698b0b4c70SToby Isaac i++; 3708b0b4c70SToby Isaac } 3718b0b4c70SToby Isaac } 3728b0b4c70SToby Isaac *adjSize = numAdj; 3738b0b4c70SToby Isaac ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 3748b0b4c70SToby Isaac } 37570034214SMatthew G. Knepley PetscFunctionReturn(0); 37670034214SMatthew G. Knepley } 37770034214SMatthew G. Knepley 37870034214SMatthew G. Knepley /*@ 37970034214SMatthew G. Knepley DMPlexGetAdjacency - Return all points adjacent to the given point 38070034214SMatthew G. Knepley 38170034214SMatthew G. Knepley Input Parameters: 38270034214SMatthew G. Knepley + dm - The DM object 38370034214SMatthew G. Knepley . p - The point 38470034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE 38570034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize 38670034214SMatthew G. Knepley 38770034214SMatthew G. Knepley Output Parameters: 38870034214SMatthew G. Knepley + adjSize - The number of adjacent points 38970034214SMatthew G. Knepley - adj - The adjacent points 39070034214SMatthew G. Knepley 39170034214SMatthew G. Knepley Level: advanced 39270034214SMatthew G. Knepley 39370034214SMatthew G. Knepley Notes: The user must PetscFree the adj array if it was not passed in. 39470034214SMatthew G. Knepley 39570034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator() 39670034214SMatthew G. Knepley @*/ 39770034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[]) 39870034214SMatthew G. Knepley { 39970034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 40070034214SMatthew G. Knepley PetscErrorCode ierr; 40170034214SMatthew G. Knepley 40270034214SMatthew G. Knepley PetscFunctionBeginHot; 40370034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 40470034214SMatthew G. Knepley PetscValidPointer(adjSize,3); 40570034214SMatthew G. Knepley PetscValidPointer(adj,4); 4065b317d89SToby Isaac ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr); 40770034214SMatthew G. Knepley PetscFunctionReturn(0); 40870034214SMatthew G. Knepley } 40908633170SToby Isaac 410b0a623aaSMatthew G. Knepley /*@ 411b0a623aaSMatthew G. Knepley DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity 412b0a623aaSMatthew G. Knepley 413b0a623aaSMatthew G. Knepley Collective on DM 414b0a623aaSMatthew G. Knepley 415b0a623aaSMatthew G. Knepley Input Parameters: 416b0a623aaSMatthew G. Knepley + dm - The DM 417b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity 418b0a623aaSMatthew G. Knepley 419b0a623aaSMatthew G. Knepley Output Parameters: 420b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL 421b0a623aaSMatthew G. Knepley - sfProcess - An SF encoding the two-sided process connectivity, or NULL 422b0a623aaSMatthew G. Knepley 423b0a623aaSMatthew G. Knepley Level: developer 424b0a623aaSMatthew G. Knepley 425b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF() 426b0a623aaSMatthew G. Knepley @*/ 427b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess) 428b0a623aaSMatthew G. Knepley { 429b0a623aaSMatthew G. Knepley const PetscSFNode *remotePoints; 430b0a623aaSMatthew G. Knepley PetscInt *localPointsNew; 431b0a623aaSMatthew G. Knepley PetscSFNode *remotePointsNew; 432b0a623aaSMatthew G. Knepley const PetscInt *nranks; 433b0a623aaSMatthew G. Knepley PetscInt *ranksNew; 434b0a623aaSMatthew G. Knepley PetscBT neighbors; 435b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, numLeaves, l, numNeighbors, n; 4369852e123SBarry Smith PetscMPIInt size, proc, rank; 437b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 438b0a623aaSMatthew G. Knepley 439b0a623aaSMatthew G. Knepley PetscFunctionBegin; 440b0a623aaSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 441b0a623aaSMatthew G. Knepley PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 442b0a623aaSMatthew G. Knepley if (processRanks) {PetscValidPointer(processRanks, 3);} 443b0a623aaSMatthew G. Knepley if (sfProcess) {PetscValidPointer(sfProcess, 4);} 4449852e123SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 445b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 446b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr); 4479852e123SBarry Smith ierr = PetscBTCreate(size, &neighbors);CHKERRQ(ierr); 4489852e123SBarry Smith ierr = PetscBTMemzero(size, neighbors);CHKERRQ(ierr); 449b0a623aaSMatthew G. Knepley /* Compute root-to-leaf process connectivity */ 450b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr); 451b0a623aaSMatthew G. Knepley ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr); 452b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 453b0a623aaSMatthew G. Knepley PetscInt ndof, noff, n; 454b0a623aaSMatthew G. Knepley 455b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr); 456b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr); 457302440fdSBarry Smith for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);} 458b0a623aaSMatthew G. Knepley } 459b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr); 460b0a623aaSMatthew G. Knepley /* Compute leaf-to-neighbor process connectivity */ 461b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr); 462b0a623aaSMatthew G. Knepley ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr); 463b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 464b0a623aaSMatthew G. Knepley PetscInt ndof, noff, n; 465b0a623aaSMatthew G. Knepley 466b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr); 467b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr); 468302440fdSBarry Smith for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);} 469b0a623aaSMatthew G. Knepley } 470b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr); 471b0a623aaSMatthew G. Knepley /* Compute leaf-to-root process connectivity */ 472b0a623aaSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);} 473b0a623aaSMatthew G. Knepley /* Calculate edges */ 474b0a623aaSMatthew G. Knepley PetscBTClear(neighbors, rank); 4759852e123SBarry Smith for(proc = 0, numNeighbors = 0; proc < size; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;} 476b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr); 477b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr); 478b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr); 4799852e123SBarry Smith for(proc = 0, n = 0; proc < size; ++proc) { 480b0a623aaSMatthew G. Knepley if (PetscBTLookup(neighbors, proc)) { 481b0a623aaSMatthew G. Knepley ranksNew[n] = proc; 482b0a623aaSMatthew G. Knepley localPointsNew[n] = proc; 483b0a623aaSMatthew G. Knepley remotePointsNew[n].index = rank; 484b0a623aaSMatthew G. Knepley remotePointsNew[n].rank = proc; 485b0a623aaSMatthew G. Knepley ++n; 486b0a623aaSMatthew G. Knepley } 487b0a623aaSMatthew G. Knepley } 488b0a623aaSMatthew G. Knepley ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr); 489b0a623aaSMatthew G. Knepley if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);} 490b0a623aaSMatthew G. Knepley else {ierr = PetscFree(ranksNew);CHKERRQ(ierr);} 491b0a623aaSMatthew G. Knepley if (sfProcess) { 492b0a623aaSMatthew G. Knepley ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr); 493b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr); 494b0a623aaSMatthew G. Knepley ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr); 4959852e123SBarry Smith ierr = PetscSFSetGraph(*sfProcess, size, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 496b0a623aaSMatthew G. Knepley } 497b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 498b0a623aaSMatthew G. Knepley } 499b0a623aaSMatthew G. Knepley 500b0a623aaSMatthew G. Knepley /*@ 501b0a623aaSMatthew G. Knepley DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF. 502b0a623aaSMatthew G. Knepley 503b0a623aaSMatthew G. Knepley Collective on DM 504b0a623aaSMatthew G. Knepley 505b0a623aaSMatthew G. Knepley Input Parameter: 506b0a623aaSMatthew G. Knepley . dm - The DM 507b0a623aaSMatthew G. Knepley 508b0a623aaSMatthew G. Knepley Output Parameters: 509b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point 510b0a623aaSMatthew G. Knepley . rootrank - The rank of each edge into the root point 511b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point 512b0a623aaSMatthew G. Knepley - leafrank - The rank of each process sharing a leaf point 513b0a623aaSMatthew G. Knepley 514b0a623aaSMatthew G. Knepley Level: developer 515b0a623aaSMatthew G. Knepley 516b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap() 517b0a623aaSMatthew G. Knepley @*/ 518b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank) 519b0a623aaSMatthew G. Knepley { 520b0a623aaSMatthew G. Knepley MPI_Comm comm; 521b0a623aaSMatthew G. Knepley PetscSF sfPoint; 522b0a623aaSMatthew G. Knepley const PetscInt *rootdegree; 523b0a623aaSMatthew G. Knepley PetscInt *myrank, *remoterank; 524b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, nedges; 525b0a623aaSMatthew G. Knepley PetscMPIInt rank; 526b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 527b0a623aaSMatthew G. Knepley 528b0a623aaSMatthew G. Knepley PetscFunctionBegin; 529b0a623aaSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 530b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 531b0a623aaSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 532b0a623aaSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 533b0a623aaSMatthew G. Knepley /* Compute number of leaves for each root */ 534b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr); 535b0a623aaSMatthew G. Knepley ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr); 536b0a623aaSMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr); 537b0a623aaSMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr); 538b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);} 539b0a623aaSMatthew G. Knepley ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 540b0a623aaSMatthew G. Knepley /* Gather rank of each leaf to root */ 541b0a623aaSMatthew G. Knepley ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr); 542b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr); 543b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(nedges, &remoterank);CHKERRQ(ierr); 544b0a623aaSMatthew G. Knepley for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank; 545b0a623aaSMatthew G. Knepley ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 546b0a623aaSMatthew G. Knepley ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 547b0a623aaSMatthew G. Knepley ierr = PetscFree(myrank);CHKERRQ(ierr); 548b0a623aaSMatthew G. Knepley ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr); 549b0a623aaSMatthew G. Knepley /* Distribute remote ranks to leaves */ 550b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr); 551b0a623aaSMatthew G. Knepley ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr); 552b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 553b0a623aaSMatthew G. Knepley } 554b0a623aaSMatthew G. Knepley 555278397a0SMatthew G. Knepley /*@C 556b0a623aaSMatthew G. Knepley DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF. 557b0a623aaSMatthew G. Knepley 558b0a623aaSMatthew G. Knepley Collective on DM 559b0a623aaSMatthew G. Knepley 560b0a623aaSMatthew G. Knepley Input Parameters: 561b0a623aaSMatthew G. Knepley + dm - The DM 56224d039d7SMichael Lange . levels - Number of overlap levels 563b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point 564b0a623aaSMatthew G. Knepley . rootrank - The rank of each edge into the root point 565b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point 566b0a623aaSMatthew G. Knepley - leafrank - The rank of each process sharing a leaf point 567b0a623aaSMatthew G. Knepley 568b0a623aaSMatthew G. Knepley Output Parameters: 569a9f1d5b2SMichael Lange + ovLabel - DMLabel containing remote overlap contributions as point/rank pairings 570b0a623aaSMatthew G. Knepley 571b0a623aaSMatthew G. Knepley Level: developer 572b0a623aaSMatthew G. Knepley 5731fd9873aSMichael Lange .seealso: DMPlexDistributeOwnership(), DMPlexDistribute() 574b0a623aaSMatthew G. Knepley @*/ 575a9f1d5b2SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel) 576b0a623aaSMatthew G. Knepley { 577e540f424SMichael Lange MPI_Comm comm; 578b0a623aaSMatthew G. Knepley DMLabel ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */ 579b0a623aaSMatthew G. Knepley PetscSF sfPoint, sfProc; 580b0a623aaSMatthew G. Knepley const PetscSFNode *remote; 581b0a623aaSMatthew G. Knepley const PetscInt *local; 5821fd9873aSMichael Lange const PetscInt *nrank, *rrank; 583b0a623aaSMatthew G. Knepley PetscInt *adj = NULL; 5841fd9873aSMichael Lange PetscInt pStart, pEnd, p, sStart, sEnd, nleaves, l; 5859852e123SBarry Smith PetscMPIInt rank, size; 58626a7d390SMatthew G. Knepley PetscBool useCone, useClosure, flg; 587b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 588b0a623aaSMatthew G. Knepley 589b0a623aaSMatthew G. Knepley PetscFunctionBegin; 590e540f424SMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 5919852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 592e540f424SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 593b0a623aaSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 594b0a623aaSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 595b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr); 596b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 597b0a623aaSMatthew G. Knepley ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr); 598b0a623aaSMatthew G. Knepley /* Handle leaves: shared with the root point */ 599b0a623aaSMatthew G. Knepley for (l = 0; l < nleaves; ++l) { 600b0a623aaSMatthew G. Knepley PetscInt adjSize = PETSC_DETERMINE, a; 601b0a623aaSMatthew G. Knepley 602b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr); 603b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);} 604b0a623aaSMatthew G. Knepley } 605b0a623aaSMatthew G. Knepley ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr); 606b0a623aaSMatthew G. Knepley ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr); 607b0a623aaSMatthew G. Knepley /* Handle roots */ 608b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 609b0a623aaSMatthew G. Knepley PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a; 610b0a623aaSMatthew G. Knepley 611b0a623aaSMatthew G. Knepley if ((p >= sStart) && (p < sEnd)) { 612b0a623aaSMatthew G. Knepley /* Some leaves share a root with other leaves on different processes */ 613b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr); 614b0a623aaSMatthew G. Knepley if (neighbors) { 615b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr); 616b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 617b0a623aaSMatthew G. Knepley for (n = 0; n < neighbors; ++n) { 618b0a623aaSMatthew G. Knepley const PetscInt remoteRank = nrank[noff+n]; 619b0a623aaSMatthew G. Knepley 620b0a623aaSMatthew G. Knepley if (remoteRank == rank) continue; 621b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 622b0a623aaSMatthew G. Knepley } 623b0a623aaSMatthew G. Knepley } 624b0a623aaSMatthew G. Knepley } 625b0a623aaSMatthew G. Knepley /* Roots are shared with leaves */ 626b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr); 627b0a623aaSMatthew G. Knepley if (!neighbors) continue; 628b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr); 629b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 630b0a623aaSMatthew G. Knepley for (n = 0; n < neighbors; ++n) { 631b0a623aaSMatthew G. Knepley const PetscInt remoteRank = rrank[noff+n]; 632b0a623aaSMatthew G. Knepley 633b0a623aaSMatthew G. Knepley if (remoteRank == rank) continue; 634b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 635b0a623aaSMatthew G. Knepley } 636b0a623aaSMatthew G. Knepley } 637b0a623aaSMatthew G. Knepley ierr = PetscFree(adj);CHKERRQ(ierr); 638b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr); 639b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr); 64024d039d7SMichael Lange /* Add additional overlap levels */ 641be200f8dSMichael Lange for (l = 1; l < levels; l++) { 642be200f8dSMichael Lange /* Propagate point donations over SF to capture remote connections */ 643be200f8dSMichael Lange ierr = DMPlexPartitionLabelPropagate(dm, ovAdjByRank);CHKERRQ(ierr); 644be200f8dSMichael Lange /* Add next level of point donations to the label */ 645be200f8dSMichael Lange ierr = DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);CHKERRQ(ierr); 646be200f8dSMichael Lange } 64726a7d390SMatthew G. Knepley /* We require the closure in the overlap */ 64826a7d390SMatthew G. Knepley ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr); 64926a7d390SMatthew G. Knepley ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr); 65026a7d390SMatthew G. Knepley if (useCone || !useClosure) { 6515abbe4feSMichael Lange ierr = DMPlexPartitionLabelClosure(dm, ovAdjByRank);CHKERRQ(ierr); 65226a7d390SMatthew G. Knepley } 653c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr); 654e540f424SMichael Lange if (flg) { 655b0a623aaSMatthew G. Knepley ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 656b0a623aaSMatthew G. Knepley } 65771a8c5fcSMichael Lange /* Make global process SF and invert sender to receiver label */ 65871a8c5fcSMichael Lange { 65971a8c5fcSMichael Lange /* Build a global process SF */ 66071a8c5fcSMichael Lange PetscSFNode *remoteProc; 6619852e123SBarry Smith ierr = PetscMalloc1(size, &remoteProc);CHKERRQ(ierr); 6629852e123SBarry Smith for (p = 0; p < size; ++p) { 66371a8c5fcSMichael Lange remoteProc[p].rank = p; 66471a8c5fcSMichael Lange remoteProc[p].index = rank; 66571a8c5fcSMichael Lange } 66671a8c5fcSMichael Lange ierr = PetscSFCreate(comm, &sfProc);CHKERRQ(ierr); 66771a8c5fcSMichael Lange ierr = PetscObjectSetName((PetscObject) sfProc, "Process SF");CHKERRQ(ierr); 6689852e123SBarry Smith ierr = PetscSFSetGraph(sfProc, size, size, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr); 66971a8c5fcSMichael Lange } 670a9f1d5b2SMichael Lange ierr = DMLabelCreate("Overlap label", ovLabel);CHKERRQ(ierr); 671a9f1d5b2SMichael Lange ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, *ovLabel);CHKERRQ(ierr); 672a9f1d5b2SMichael Lange /* Add owned points, except for shared local points */ 673a9f1d5b2SMichael Lange for (p = pStart; p < pEnd; ++p) {ierr = DMLabelSetValue(*ovLabel, p, rank);CHKERRQ(ierr);} 674e540f424SMichael Lange for (l = 0; l < nleaves; ++l) { 675a9f1d5b2SMichael Lange ierr = DMLabelClearValue(*ovLabel, local[l], rank);CHKERRQ(ierr); 676a9f1d5b2SMichael Lange ierr = DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr); 677e540f424SMichael Lange } 678e540f424SMichael Lange /* Clean up */ 6791fd9873aSMichael Lange ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr); 680b0a623aaSMatthew G. Knepley ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr); 681b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 682b0a623aaSMatthew G. Knepley } 68370034214SMatthew G. Knepley 68424cc2ca5SMatthew G. Knepley /*@C 68524cc2ca5SMatthew G. Knepley DMPlexCreateOverlapMigrationSF - Create an SF describing the new mesh distribution to make the overlap described by the input SF 68624cc2ca5SMatthew G. Knepley 68724cc2ca5SMatthew G. Knepley Collective on DM 68824cc2ca5SMatthew G. Knepley 68924cc2ca5SMatthew G. Knepley Input Parameters: 69024cc2ca5SMatthew G. Knepley + dm - The DM 69124cc2ca5SMatthew G. Knepley - overlapSF - The SF mapping ghost points in overlap to owner points on other processes 69224cc2ca5SMatthew G. Knepley 69324cc2ca5SMatthew G. Knepley Output Parameters: 69424cc2ca5SMatthew G. Knepley + migrationSF - An SF that maps original points in old locations to points in new locations 69524cc2ca5SMatthew G. Knepley 69624cc2ca5SMatthew G. Knepley Level: developer 69724cc2ca5SMatthew G. Knepley 69824cc2ca5SMatthew G. Knepley .seealso: DMPlexCreateOverlap(), DMPlexDistribute() 69924cc2ca5SMatthew G. Knepley @*/ 70046f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF) 70146f9b1c3SMichael Lange { 70246f9b1c3SMichael Lange MPI_Comm comm; 7039852e123SBarry Smith PetscMPIInt rank, size; 70446f9b1c3SMichael Lange PetscInt d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints; 70546f9b1c3SMichael Lange PetscInt *pointDepths, *remoteDepths, *ilocal; 70646f9b1c3SMichael Lange PetscInt *depthRecv, *depthShift, *depthIdx; 70746f9b1c3SMichael Lange PetscSFNode *iremote; 70846f9b1c3SMichael Lange PetscSF pointSF; 70946f9b1c3SMichael Lange const PetscInt *sharedLocal; 71046f9b1c3SMichael Lange const PetscSFNode *overlapRemote, *sharedRemote; 71146f9b1c3SMichael Lange PetscErrorCode ierr; 71246f9b1c3SMichael Lange 71346f9b1c3SMichael Lange PetscFunctionBegin; 71446f9b1c3SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 71546f9b1c3SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 71646f9b1c3SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 7179852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 71846f9b1c3SMichael Lange ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 71946f9b1c3SMichael Lange 72046f9b1c3SMichael Lange /* Before building the migration SF we need to know the new stratum offsets */ 72146f9b1c3SMichael Lange ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr); 72246f9b1c3SMichael Lange ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr); 72346f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 72446f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 72546f9b1c3SMichael Lange for (p=pStart; p<pEnd; p++) pointDepths[p] = d; 72646f9b1c3SMichael Lange } 72746f9b1c3SMichael Lange for (p=0; p<nleaves; p++) remoteDepths[p] = -1; 72846f9b1c3SMichael Lange ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 72946f9b1c3SMichael Lange ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 73046f9b1c3SMichael Lange 73146f9b1c3SMichael Lange /* Count recevied points in each stratum and compute the internal strata shift */ 73246f9b1c3SMichael Lange ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr); 73346f9b1c3SMichael Lange for (d=0; d<dim+1; d++) depthRecv[d]=0; 73446f9b1c3SMichael Lange for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++; 73546f9b1c3SMichael Lange depthShift[dim] = 0; 73646f9b1c3SMichael Lange for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim]; 73746f9b1c3SMichael Lange for (d=1; d<dim; d++) depthShift[d] += depthRecv[0]; 73846f9b1c3SMichael Lange for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1]; 73946f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 74046f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 74146f9b1c3SMichael Lange depthIdx[d] = pStart + depthShift[d]; 74246f9b1c3SMichael Lange } 74346f9b1c3SMichael Lange 74446f9b1c3SMichael Lange /* Form the overlap SF build an SF that describes the full overlap migration SF */ 74546f9b1c3SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 74646f9b1c3SMichael Lange newLeaves = pEnd - pStart + nleaves; 74709b7985cSMatthew G. Knepley ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr); 74809b7985cSMatthew G. Knepley ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr); 74946f9b1c3SMichael Lange /* First map local points to themselves */ 75046f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 75146f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 75246f9b1c3SMichael Lange for (p=pStart; p<pEnd; p++) { 75346f9b1c3SMichael Lange point = p + depthShift[d]; 75446f9b1c3SMichael Lange ilocal[point] = point; 75546f9b1c3SMichael Lange iremote[point].index = p; 75646f9b1c3SMichael Lange iremote[point].rank = rank; 75746f9b1c3SMichael Lange depthIdx[d]++; 75846f9b1c3SMichael Lange } 75946f9b1c3SMichael Lange } 76046f9b1c3SMichael Lange 76146f9b1c3SMichael Lange /* Add in the remote roots for currently shared points */ 76246f9b1c3SMichael Lange ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr); 76346f9b1c3SMichael Lange ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr); 76446f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 76546f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 76646f9b1c3SMichael Lange for (p=0; p<numSharedPoints; p++) { 76746f9b1c3SMichael Lange if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) { 76846f9b1c3SMichael Lange point = sharedLocal[p] + depthShift[d]; 76946f9b1c3SMichael Lange iremote[point].index = sharedRemote[p].index; 77046f9b1c3SMichael Lange iremote[point].rank = sharedRemote[p].rank; 77146f9b1c3SMichael Lange } 77246f9b1c3SMichael Lange } 77346f9b1c3SMichael Lange } 77446f9b1c3SMichael Lange 77546f9b1c3SMichael Lange /* Now add the incoming overlap points */ 77646f9b1c3SMichael Lange for (p=0; p<nleaves; p++) { 77746f9b1c3SMichael Lange point = depthIdx[remoteDepths[p]]; 77846f9b1c3SMichael Lange ilocal[point] = point; 77946f9b1c3SMichael Lange iremote[point].index = overlapRemote[p].index; 78046f9b1c3SMichael Lange iremote[point].rank = overlapRemote[p].rank; 78146f9b1c3SMichael Lange depthIdx[remoteDepths[p]]++; 78246f9b1c3SMichael Lange } 78315fff7beSMatthew G. Knepley ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr); 78446f9b1c3SMichael Lange 78546f9b1c3SMichael Lange ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr); 78646f9b1c3SMichael Lange ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr); 78746f9b1c3SMichael Lange ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr); 78846f9b1c3SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 78946f9b1c3SMichael Lange ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 79046f9b1c3SMichael Lange 79146f9b1c3SMichael Lange ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr); 79270034214SMatthew G. Knepley PetscFunctionReturn(0); 79370034214SMatthew G. Knepley } 79470034214SMatthew G. Knepley 795a9f1d5b2SMichael Lange /*@ 796f0e73a3dSToby Isaac DMPlexStratifyMigrationSF - Rearrange the leaves of a migration sf for stratification. 797a9f1d5b2SMichael Lange 798a9f1d5b2SMichael Lange Input Parameter: 799a9f1d5b2SMichael Lange + dm - The DM 800a9f1d5b2SMichael Lange - sf - A star forest with non-ordered leaves, usually defining a DM point migration 801a9f1d5b2SMichael Lange 802a9f1d5b2SMichael Lange Output Parameter: 803a9f1d5b2SMichael Lange . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified 804a9f1d5b2SMichael Lange 805a9f1d5b2SMichael Lange Level: developer 806a9f1d5b2SMichael Lange 807a9f1d5b2SMichael Lange .seealso: DMPlexPartitionLabelCreateSF(), DMPlexDistribute(), DMPlexDistributeOverlap() 808a9f1d5b2SMichael Lange @*/ 809a9f1d5b2SMichael Lange PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF) 810a9f1d5b2SMichael Lange { 811a9f1d5b2SMichael Lange MPI_Comm comm; 8129852e123SBarry Smith PetscMPIInt rank, size; 8137fab53ddSMatthew G. Knepley PetscInt d, ldepth, depth, p, pStart, pEnd, nroots, nleaves; 814a9f1d5b2SMichael Lange PetscInt *pointDepths, *remoteDepths, *ilocal; 815a9f1d5b2SMichael Lange PetscInt *depthRecv, *depthShift, *depthIdx; 816f0e73a3dSToby Isaac PetscInt hybEnd[4]; 817a9f1d5b2SMichael Lange const PetscSFNode *iremote; 818a9f1d5b2SMichael Lange PetscErrorCode ierr; 819a9f1d5b2SMichael Lange 820a9f1d5b2SMichael Lange PetscFunctionBegin; 821a9f1d5b2SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 822a9f1d5b2SMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 823a9f1d5b2SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 8249852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 8257fab53ddSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &ldepth);CHKERRQ(ierr); 826b2566f29SBarry Smith ierr = MPIU_Allreduce(&ldepth, &depth, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 8277fab53ddSMatthew G. Knepley if ((ldepth >= 0) && (depth != ldepth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", ldepth, depth); 828a9f1d5b2SMichael Lange 829a9f1d5b2SMichael Lange /* Before building the migration SF we need to know the new stratum offsets */ 830a9f1d5b2SMichael Lange ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote);CHKERRQ(ierr); 831a9f1d5b2SMichael Lange ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr); 832f0e73a3dSToby Isaac ierr = DMPlexGetHybridBounds(dm,&hybEnd[depth],&hybEnd[depth-1],&hybEnd[1],&hybEnd[0]);CHKERRQ(ierr); 8337fab53ddSMatthew G. Knepley for (d = 0; d < depth+1; ++d) { 834a9f1d5b2SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 835f0e73a3dSToby Isaac for (p = pStart; p < pEnd; ++p) { 836f0e73a3dSToby Isaac if (hybEnd[d] >= 0 && p >= hybEnd[d]) { /* put in a separate value for hybrid points */ 837f0e73a3dSToby Isaac pointDepths[p] = 2 * d; 838f0e73a3dSToby Isaac } else { 839f0e73a3dSToby Isaac pointDepths[p] = 2 * d + 1; 840f0e73a3dSToby Isaac } 841f0e73a3dSToby Isaac } 842a9f1d5b2SMichael Lange } 8437fab53ddSMatthew G. Knepley for (p = 0; p < nleaves; ++p) remoteDepths[p] = -1; 844a9f1d5b2SMichael Lange ierr = PetscSFBcastBegin(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 845a9f1d5b2SMichael Lange ierr = PetscSFBcastEnd(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 846a9f1d5b2SMichael Lange /* Count recevied points in each stratum and compute the internal strata shift */ 847f0e73a3dSToby Isaac ierr = PetscMalloc3(2*(depth+1), &depthRecv, 2*(depth+1), &depthShift, 2*(depth+1), &depthIdx);CHKERRQ(ierr); 848f0e73a3dSToby Isaac for (d = 0; d < 2*(depth+1); ++d) depthRecv[d] = 0; 8497fab53ddSMatthew G. Knepley for (p = 0; p < nleaves; ++p) depthRecv[remoteDepths[p]]++; 850f0e73a3dSToby Isaac depthShift[2*depth+1] = 0; 851f0e73a3dSToby Isaac for (d = 0; d < 2*depth+1; ++d) depthShift[d] = depthRecv[2 * depth + 1]; 852f0e73a3dSToby Isaac for (d = 0; d < 2*depth; ++d) depthShift[d] += depthRecv[2 * depth]; 853f0e73a3dSToby Isaac depthShift[0] += depthRecv[1]; 854f0e73a3dSToby Isaac for (d = 2; d < 2*depth; ++d) depthShift[d] += depthRecv[1]; 855f0e73a3dSToby Isaac for (d = 2; d < 2*depth; ++d) depthShift[d] += depthRecv[0]; 856f0e73a3dSToby Isaac for (d = 2 * depth-1; d > 2; --d) { 857f0e73a3dSToby Isaac PetscInt e; 858f0e73a3dSToby Isaac 859f0e73a3dSToby Isaac for (e = d -1; e > 1; --e) depthShift[e] += depthRecv[d]; 860f0e73a3dSToby Isaac } 861f0e73a3dSToby Isaac for (d = 0; d < 2*(depth+1); ++d) {depthIdx[d] = 0;} 862a9f1d5b2SMichael Lange /* Derive a new local permutation based on stratified indices */ 863a9f1d5b2SMichael Lange ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 8647fab53ddSMatthew G. Knepley for (p = 0; p < nleaves; ++p) { 8657fab53ddSMatthew G. Knepley const PetscInt dep = remoteDepths[p]; 8667fab53ddSMatthew G. Knepley 8677fab53ddSMatthew G. Knepley ilocal[p] = depthShift[dep] + depthIdx[dep]; 8687fab53ddSMatthew G. Knepley depthIdx[dep]++; 869a9f1d5b2SMichael Lange } 870a9f1d5b2SMichael Lange ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr); 871a9f1d5b2SMichael Lange ierr = PetscObjectSetName((PetscObject) *migrationSF, "Migration SF");CHKERRQ(ierr); 872a9f1d5b2SMichael Lange ierr = PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_COPY_VALUES);CHKERRQ(ierr); 873a9f1d5b2SMichael Lange ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr); 874a9f1d5b2SMichael Lange ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr); 875a9f1d5b2SMichael Lange PetscFunctionReturn(0); 876a9f1d5b2SMichael Lange } 877a9f1d5b2SMichael Lange 87870034214SMatthew G. Knepley /*@ 87970034214SMatthew G. Knepley DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 88070034214SMatthew G. Knepley 88170034214SMatthew G. Knepley Collective on DM 88270034214SMatthew G. Knepley 88370034214SMatthew G. Knepley Input Parameters: 88470034214SMatthew G. Knepley + dm - The DMPlex object 88570034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 88670034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 88770034214SMatthew G. Knepley - originalVec - The existing data 88870034214SMatthew G. Knepley 88970034214SMatthew G. Knepley Output Parameters: 89070034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 89170034214SMatthew G. Knepley - newVec - The new data 89270034214SMatthew G. Knepley 89370034214SMatthew G. Knepley Level: developer 89470034214SMatthew G. Knepley 89570034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData() 89670034214SMatthew G. Knepley @*/ 89770034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec) 89870034214SMatthew G. Knepley { 89970034214SMatthew G. Knepley PetscSF fieldSF; 90070034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 90170034214SMatthew G. Knepley PetscScalar *originalValues, *newValues; 90270034214SMatthew G. Knepley PetscErrorCode ierr; 90370034214SMatthew G. Knepley 90470034214SMatthew G. Knepley PetscFunctionBegin; 90570034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 90670034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 90770034214SMatthew G. Knepley 90870034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 90970034214SMatthew G. Knepley ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr); 91070034214SMatthew G. Knepley ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr); 91170034214SMatthew G. Knepley 91270034214SMatthew G. Knepley ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr); 91370034214SMatthew G. Knepley ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr); 91470034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 9158a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 91670034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 91770034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 91870034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 91970034214SMatthew G. Knepley ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr); 92070034214SMatthew G. Knepley ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr); 92170034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 92270034214SMatthew G. Knepley PetscFunctionReturn(0); 92370034214SMatthew G. Knepley } 92470034214SMatthew G. Knepley 92570034214SMatthew G. Knepley /*@ 92670034214SMatthew G. Knepley DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 92770034214SMatthew G. Knepley 92870034214SMatthew G. Knepley Collective on DM 92970034214SMatthew G. Knepley 93070034214SMatthew G. Knepley Input Parameters: 93170034214SMatthew G. Knepley + dm - The DMPlex object 93270034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 93370034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 93470034214SMatthew G. Knepley - originalIS - The existing data 93570034214SMatthew G. Knepley 93670034214SMatthew G. Knepley Output Parameters: 93770034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 93870034214SMatthew G. Knepley - newIS - The new data 93970034214SMatthew G. Knepley 94070034214SMatthew G. Knepley Level: developer 94170034214SMatthew G. Knepley 94270034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData() 94370034214SMatthew G. Knepley @*/ 94470034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS) 94570034214SMatthew G. Knepley { 94670034214SMatthew G. Knepley PetscSF fieldSF; 94770034214SMatthew G. Knepley PetscInt *newValues, *remoteOffsets, fieldSize; 94870034214SMatthew G. Knepley const PetscInt *originalValues; 94970034214SMatthew G. Knepley PetscErrorCode ierr; 95070034214SMatthew G. Knepley 95170034214SMatthew G. Knepley PetscFunctionBegin; 95270034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 95370034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 95470034214SMatthew G. Knepley 95570034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 956854ce69bSBarry Smith ierr = PetscMalloc1(fieldSize, &newValues);CHKERRQ(ierr); 95770034214SMatthew G. Knepley 95870034214SMatthew G. Knepley ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr); 95970034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 9608a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 961aaf8c182SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr); 962aaf8c182SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr); 96370034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 96470034214SMatthew G. Knepley ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr); 96570034214SMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr); 96670034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 96770034214SMatthew G. Knepley PetscFunctionReturn(0); 96870034214SMatthew G. Knepley } 96970034214SMatthew G. Knepley 97070034214SMatthew G. Knepley /*@ 97170034214SMatthew G. Knepley DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 97270034214SMatthew G. Knepley 97370034214SMatthew G. Knepley Collective on DM 97470034214SMatthew G. Knepley 97570034214SMatthew G. Knepley Input Parameters: 97670034214SMatthew G. Knepley + dm - The DMPlex object 97770034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 97870034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 97970034214SMatthew G. Knepley . datatype - The type of data 98070034214SMatthew G. Knepley - originalData - The existing data 98170034214SMatthew G. Knepley 98270034214SMatthew G. Knepley Output Parameters: 98370034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout 98470034214SMatthew G. Knepley - newData - The new data 98570034214SMatthew G. Knepley 98670034214SMatthew G. Knepley Level: developer 98770034214SMatthew G. Knepley 98870034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField() 98970034214SMatthew G. Knepley @*/ 99070034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData) 99170034214SMatthew G. Knepley { 99270034214SMatthew G. Knepley PetscSF fieldSF; 99370034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 99470034214SMatthew G. Knepley PetscMPIInt dataSize; 99570034214SMatthew G. Knepley PetscErrorCode ierr; 99670034214SMatthew G. Knepley 99770034214SMatthew G. Knepley PetscFunctionBegin; 99870034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 99970034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 100070034214SMatthew G. Knepley 100170034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 100270034214SMatthew G. Knepley ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr); 100370034214SMatthew G. Knepley ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr); 100470034214SMatthew G. Knepley 100570034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 10068a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 100770034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 100870034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 100970034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 101070034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 101170034214SMatthew G. Knepley PetscFunctionReturn(0); 101270034214SMatthew G. Knepley } 101370034214SMatthew G. Knepley 101424cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel) 1015cc71bff1SMichael Lange { 1016f5bf2dbfSMichael Lange DM_Plex *mesh = (DM_Plex*) dm->data; 1017cc71bff1SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 1018cc71bff1SMichael Lange MPI_Comm comm; 1019cc71bff1SMichael Lange PetscSF coneSF; 1020cc71bff1SMichael Lange PetscSection originalConeSection, newConeSection; 1021ac04eaf7SToby Isaac PetscInt *remoteOffsets, *cones, *globCones, *newCones, newConesSize; 1022cc71bff1SMichael Lange PetscBool flg; 1023cc71bff1SMichael Lange PetscErrorCode ierr; 1024cc71bff1SMichael Lange 1025cc71bff1SMichael Lange PetscFunctionBegin; 1026cc71bff1SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1027cc71bff1SMichael Lange PetscValidPointer(dmParallel,4); 1028cc71bff1SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 1029cc71bff1SMichael Lange 1030cc71bff1SMichael Lange /* Distribute cone section */ 1031cc71bff1SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1032cc71bff1SMichael Lange ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); 1033cc71bff1SMichael Lange ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr); 1034cc71bff1SMichael Lange ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); 1035cc71bff1SMichael Lange ierr = DMSetUp(dmParallel);CHKERRQ(ierr); 1036cc71bff1SMichael Lange { 1037cc71bff1SMichael Lange PetscInt pStart, pEnd, p; 1038cc71bff1SMichael Lange 1039cc71bff1SMichael Lange ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr); 1040cc71bff1SMichael Lange for (p = pStart; p < pEnd; ++p) { 1041cc71bff1SMichael Lange PetscInt coneSize; 1042cc71bff1SMichael Lange ierr = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr); 1043cc71bff1SMichael Lange pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize); 1044cc71bff1SMichael Lange } 1045cc71bff1SMichael Lange } 1046cc71bff1SMichael Lange /* Communicate and renumber cones */ 1047cc71bff1SMichael Lange ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); 10488a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1049cc71bff1SMichael Lange ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 1050ac04eaf7SToby Isaac if (original) { 1051ac04eaf7SToby Isaac PetscInt numCones; 1052ac04eaf7SToby Isaac 1053ac04eaf7SToby Isaac ierr = PetscSectionGetStorageSize(originalConeSection,&numCones);CHKERRQ(ierr); ierr = PetscMalloc1(numCones,&globCones);CHKERRQ(ierr); 1054ac04eaf7SToby Isaac ierr = ISLocalToGlobalMappingApplyBlock(original, numCones, cones, globCones);CHKERRQ(ierr); 1055ac04eaf7SToby Isaac } 1056ac04eaf7SToby Isaac else { 1057ac04eaf7SToby Isaac globCones = cones; 1058ac04eaf7SToby Isaac } 1059cc71bff1SMichael Lange ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr); 1060ac04eaf7SToby Isaac ierr = PetscSFBcastBegin(coneSF, MPIU_INT, globCones, newCones);CHKERRQ(ierr); 1061ac04eaf7SToby Isaac ierr = PetscSFBcastEnd(coneSF, MPIU_INT, globCones, newCones);CHKERRQ(ierr); 1062ac04eaf7SToby Isaac if (original) { 1063ac04eaf7SToby Isaac ierr = PetscFree(globCones);CHKERRQ(ierr); 1064ac04eaf7SToby Isaac } 1065cc71bff1SMichael Lange ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); 1066cc71bff1SMichael Lange ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr); 10673533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG 10683533c52bSMatthew G. Knepley { 10693533c52bSMatthew G. Knepley PetscInt p; 10703533c52bSMatthew G. Knepley PetscBool valid = PETSC_TRUE; 10713533c52bSMatthew G. Knepley for (p = 0; p < newConesSize; ++p) { 10723533c52bSMatthew G. Knepley if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);} 10733533c52bSMatthew G. Knepley } 10743533c52bSMatthew G. Knepley if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map"); 10753533c52bSMatthew G. Knepley } 10763533c52bSMatthew G. Knepley #endif 1077c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr); 1078cc71bff1SMichael Lange if (flg) { 1079cc71bff1SMichael Lange ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); 1080cc71bff1SMichael Lange ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1081cc71bff1SMichael Lange ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); 1082cc71bff1SMichael Lange ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1083cc71bff1SMichael Lange ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr); 1084cc71bff1SMichael Lange } 1085cc71bff1SMichael Lange ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr); 1086cc71bff1SMichael Lange ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr); 1087cc71bff1SMichael Lange ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 1088cc71bff1SMichael Lange ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 1089cc71bff1SMichael Lange ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); 1090cc71bff1SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 1091cc71bff1SMichael Lange /* Create supports and stratify sieve */ 1092cc71bff1SMichael Lange { 1093cc71bff1SMichael Lange PetscInt pStart, pEnd; 1094cc71bff1SMichael Lange 1095cc71bff1SMichael Lange ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 1096cc71bff1SMichael Lange ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr); 1097cc71bff1SMichael Lange } 1098cc71bff1SMichael Lange ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr); 1099cc71bff1SMichael Lange ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr); 1100f5bf2dbfSMichael Lange pmesh->useCone = mesh->useCone; 1101f5bf2dbfSMichael Lange pmesh->useClosure = mesh->useClosure; 1102cc71bff1SMichael Lange PetscFunctionReturn(0); 1103cc71bff1SMichael Lange } 1104cc71bff1SMichael Lange 110524cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel) 11060df0e737SMichael Lange { 11070df0e737SMichael Lange MPI_Comm comm; 11080df0e737SMichael Lange PetscSection originalCoordSection, newCoordSection; 11090df0e737SMichael Lange Vec originalCoordinates, newCoordinates; 11100df0e737SMichael Lange PetscInt bs; 111190b157c4SStefano Zampini PetscBool isper; 11120df0e737SMichael Lange const char *name; 11130df0e737SMichael Lange const PetscReal *maxCell, *L; 11145dc8c3f7SMatthew G. Knepley const DMBoundaryType *bd; 11150df0e737SMichael Lange PetscErrorCode ierr; 11160df0e737SMichael Lange 11170df0e737SMichael Lange PetscFunctionBegin; 11180df0e737SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11190df0e737SMichael Lange PetscValidPointer(dmParallel, 3); 11200df0e737SMichael Lange 11210df0e737SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 11220df0e737SMichael Lange ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); 11230df0e737SMichael Lange ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr); 11240df0e737SMichael Lange ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr); 11250df0e737SMichael Lange if (originalCoordinates) { 11268b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &newCoordinates);CHKERRQ(ierr); 11270df0e737SMichael Lange ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr); 11280df0e737SMichael Lange ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr); 11290df0e737SMichael Lange 11300df0e737SMichael Lange ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr); 11310df0e737SMichael Lange ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr); 11320df0e737SMichael Lange ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr); 11330df0e737SMichael Lange ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr); 11340df0e737SMichael Lange ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 11350df0e737SMichael Lange } 113690b157c4SStefano Zampini ierr = DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);CHKERRQ(ierr); 113790b157c4SStefano Zampini ierr = DMSetPeriodicity(dmParallel, isper, maxCell, L, bd);CHKERRQ(ierr); 11380df0e737SMichael Lange PetscFunctionReturn(0); 11390df0e737SMichael Lange } 11400df0e737SMichael Lange 1141d995df53SMatthew G. Knepley /* Here we are assuming that process 0 always has everything */ 114224cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel) 11430df0e737SMichael Lange { 1144df0420ecSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 11450df0e737SMichael Lange MPI_Comm comm; 11467980c9d4SMatthew G. Knepley DMLabel depthLabel; 11470df0e737SMichael Lange PetscMPIInt rank; 11487980c9d4SMatthew G. Knepley PetscInt depth, d, numLabels, numLocalLabels, l; 1149df0420ecSMatthew G. Knepley PetscBool hasLabels = PETSC_FALSE, lsendDepth, sendDepth; 1150df0420ecSMatthew G. Knepley PetscObjectState depthState = -1; 11510df0e737SMichael Lange PetscErrorCode ierr; 11520df0e737SMichael Lange 11530df0e737SMichael Lange PetscFunctionBegin; 11540df0e737SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11552eb1fa14SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 11560df0e737SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 11570df0e737SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 11580df0e737SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 11590df0e737SMichael Lange 1160df0420ecSMatthew G. Knepley /* If the user has changed the depth label, communicate it instead */ 11617980c9d4SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 11627980c9d4SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 11637980c9d4SMatthew G. Knepley if (depthLabel) {ierr = DMLabelGetState(depthLabel, &depthState);CHKERRQ(ierr);} 1164df0420ecSMatthew G. Knepley lsendDepth = mesh->depthState != depthState ? PETSC_TRUE : PETSC_FALSE; 1165b2566f29SBarry Smith ierr = MPIU_Allreduce(&lsendDepth, &sendDepth, 1, MPIU_BOOL, MPI_LOR, comm);CHKERRQ(ierr); 1166df0420ecSMatthew G. Knepley if (sendDepth) { 11677980c9d4SMatthew G. Knepley ierr = DMRemoveLabel(dmParallel, "depth", &depthLabel);CHKERRQ(ierr); 11687980c9d4SMatthew G. Knepley ierr = DMLabelDestroy(&depthLabel);CHKERRQ(ierr); 1169df0420ecSMatthew G. Knepley } 1170d995df53SMatthew G. Knepley /* Everyone must have either the same number of labels, or none */ 1171c58f1c22SToby Isaac ierr = DMGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr); 1172d995df53SMatthew G. Knepley numLabels = numLocalLabels; 11730df0e737SMichael Lange ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1174627847f0SMatthew G. Knepley if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE; 1175bdd2d751SMichael Lange for (l = numLabels-1; l >= 0; --l) { 1176bdd2d751SMichael Lange DMLabel label = NULL, labelNew = NULL; 11770df0e737SMichael Lange PetscBool isdepth; 11780df0e737SMichael Lange 1179d995df53SMatthew G. Knepley if (hasLabels) { 1180c58f1c22SToby Isaac ierr = DMGetLabelByNum(dm, l, &label);CHKERRQ(ierr); 11810df0e737SMichael Lange /* Skip "depth" because it is recreated */ 1182bdd2d751SMichael Lange ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr); 1183bdd2d751SMichael Lange } 11840df0e737SMichael Lange ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 1185df0420ecSMatthew G. Knepley if (isdepth && !sendDepth) continue; 1186bdd2d751SMichael Lange ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr); 11877980c9d4SMatthew G. Knepley if (isdepth) { 11887980c9d4SMatthew G. Knepley /* Put in any missing strata which can occur if users are managing the depth label themselves */ 11897980c9d4SMatthew G. Knepley PetscInt gdepth; 11907980c9d4SMatthew G. Knepley 11917980c9d4SMatthew G. Knepley ierr = MPIU_Allreduce(&depth, &gdepth, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 11927980c9d4SMatthew G. Knepley if ((depth >= 0) && (gdepth != depth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", depth, gdepth); 11937980c9d4SMatthew G. Knepley for (d = 0; d <= gdepth; ++d) { 11947980c9d4SMatthew G. Knepley PetscBool has; 11957980c9d4SMatthew G. Knepley 11967980c9d4SMatthew G. Knepley ierr = DMLabelHasStratum(labelNew, d, &has);CHKERRQ(ierr); 11977980c9d4SMatthew G. Knepley if (!has) ierr = DMLabelAddStratum(labelNew, d);CHKERRQ(ierr); 11987980c9d4SMatthew G. Knepley } 11997980c9d4SMatthew G. Knepley } 1200c58f1c22SToby Isaac ierr = DMAddLabel(dmParallel, labelNew);CHKERRQ(ierr); 12010df0e737SMichael Lange } 12020df0e737SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 12030df0e737SMichael Lange PetscFunctionReturn(0); 12040df0e737SMichael Lange } 12050df0e737SMichael Lange 120624cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 12079734c634SMichael Lange { 12089734c634SMichael Lange DM_Plex *mesh = (DM_Plex*) dm->data; 12099734c634SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 1210f0e73a3dSToby Isaac PetscBool *isHybrid, *isHybridParallel; 1211f0e73a3dSToby Isaac PetscInt dim, depth, d; 1212f0e73a3dSToby Isaac PetscInt pStart, pEnd, pStartP, pEndP; 12139734c634SMichael Lange PetscErrorCode ierr; 12149734c634SMichael Lange 12159734c634SMichael Lange PetscFunctionBegin; 12169734c634SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12179734c634SMichael Lange PetscValidPointer(dmParallel, 4); 12189734c634SMichael Lange 12199734c634SMichael Lange ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 12209734c634SMichael Lange ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1221f0e73a3dSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 1222f0e73a3dSToby Isaac ierr = DMPlexGetChart(dmParallel,&pStartP,&pEndP);CHKERRQ(ierr); 1223f0e73a3dSToby Isaac ierr = PetscCalloc2(pEnd-pStart,&isHybrid,pEndP-pStartP,&isHybridParallel);CHKERRQ(ierr); 1224f0e73a3dSToby Isaac for (d = 0; d <= depth; d++) { 1225f0e73a3dSToby Isaac PetscInt hybridMax = (depth == 1 && d == 1) ? mesh->hybridPointMax[dim] : mesh->hybridPointMax[d]; 12269734c634SMichael Lange 1227f0e73a3dSToby Isaac if (hybridMax >= 0) { 1228f0e73a3dSToby Isaac PetscInt sStart, sEnd, p; 12299734c634SMichael Lange 1230f0e73a3dSToby Isaac ierr = DMPlexGetDepthStratum(dm,d,&sStart,&sEnd);CHKERRQ(ierr); 1231f0e73a3dSToby Isaac for (p = hybridMax; p < sEnd; p++) isHybrid[p-pStart] = PETSC_TRUE; 12329734c634SMichael Lange } 12339734c634SMichael Lange } 1234f0e73a3dSToby Isaac ierr = PetscSFBcastBegin(migrationSF,MPIU_BOOL,isHybrid,isHybridParallel);CHKERRQ(ierr); 1235f0e73a3dSToby Isaac ierr = PetscSFBcastEnd(migrationSF,MPIU_BOOL,isHybrid,isHybridParallel);CHKERRQ(ierr); 1236f0e73a3dSToby Isaac for (d = 0; d <= dim; d++) pmesh->hybridPointMax[d] = -1; 1237f0e73a3dSToby Isaac for (d = 0; d <= depth; d++) { 1238f0e73a3dSToby Isaac PetscInt sStart, sEnd, p, dd; 1239f0e73a3dSToby Isaac 1240f0e73a3dSToby Isaac ierr = DMPlexGetDepthStratum(dmParallel,d,&sStart,&sEnd);CHKERRQ(ierr); 1241f0e73a3dSToby Isaac dd = (depth == 1 && d == 1) ? dim : d; 1242f0e73a3dSToby Isaac for (p = sStart; p < sEnd; p++) { 1243f0e73a3dSToby Isaac if (isHybridParallel[p-pStartP]) { 1244f0e73a3dSToby Isaac pmesh->hybridPointMax[dd] = p; 1245f0e73a3dSToby Isaac break; 1246f0e73a3dSToby Isaac } 1247f0e73a3dSToby Isaac } 1248f0e73a3dSToby Isaac } 1249f0e73a3dSToby Isaac ierr = PetscFree2(isHybrid,isHybridParallel);CHKERRQ(ierr); 12509734c634SMichael Lange PetscFunctionReturn(0); 12519734c634SMichael Lange } 12520df0e737SMichael Lange 125324cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel) 1254a6f36705SMichael Lange { 125515078cd4SMichael Lange DM_Plex *mesh = (DM_Plex*) dm->data; 125615078cd4SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 1257a6f36705SMichael Lange MPI_Comm comm; 1258a6f36705SMichael Lange DM refTree; 1259a6f36705SMichael Lange PetscSection origParentSection, newParentSection; 1260a6f36705SMichael Lange PetscInt *origParents, *origChildIDs; 1261a6f36705SMichael Lange PetscBool flg; 1262a6f36705SMichael Lange PetscErrorCode ierr; 1263a6f36705SMichael Lange 1264a6f36705SMichael Lange PetscFunctionBegin; 1265a6f36705SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1266a6f36705SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 4); 1267a6f36705SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1268a6f36705SMichael Lange 1269a6f36705SMichael Lange /* Set up tree */ 1270a6f36705SMichael Lange ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 1271a6f36705SMichael Lange ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr); 1272a6f36705SMichael Lange ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr); 1273a6f36705SMichael Lange if (origParentSection) { 1274a6f36705SMichael Lange PetscInt pStart, pEnd; 127508633170SToby Isaac PetscInt *newParents, *newChildIDs, *globParents; 1276a6f36705SMichael Lange PetscInt *remoteOffsetsParents, newParentSize; 1277a6f36705SMichael Lange PetscSF parentSF; 1278a6f36705SMichael Lange 1279a6f36705SMichael Lange ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 1280a6f36705SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr); 1281a6f36705SMichael Lange ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr); 1282a6f36705SMichael Lange ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr); 1283a6f36705SMichael Lange ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr); 12848a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsetsParents);CHKERRQ(ierr); 1285a6f36705SMichael Lange ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr); 1286a6f36705SMichael Lange ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr); 128708633170SToby Isaac if (original) { 128808633170SToby Isaac PetscInt numParents; 128908633170SToby Isaac 129008633170SToby Isaac ierr = PetscSectionGetStorageSize(origParentSection,&numParents);CHKERRQ(ierr); 129108633170SToby Isaac ierr = PetscMalloc1(numParents,&globParents);CHKERRQ(ierr); 129208633170SToby Isaac ierr = ISLocalToGlobalMappingApplyBlock(original, numParents, origParents, globParents);CHKERRQ(ierr); 129308633170SToby Isaac } 129408633170SToby Isaac else { 129508633170SToby Isaac globParents = origParents; 129608633170SToby Isaac } 129708633170SToby Isaac ierr = PetscSFBcastBegin(parentSF, MPIU_INT, globParents, newParents);CHKERRQ(ierr); 129808633170SToby Isaac ierr = PetscSFBcastEnd(parentSF, MPIU_INT, globParents, newParents);CHKERRQ(ierr); 129908633170SToby Isaac if (original) { 130008633170SToby Isaac ierr = PetscFree(globParents);CHKERRQ(ierr); 130108633170SToby Isaac } 1302a6f36705SMichael Lange ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr); 1303a6f36705SMichael Lange ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr); 1304a6f36705SMichael Lange ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr); 13054a54e071SToby Isaac #if PETSC_USE_DEBUG 13064a54e071SToby Isaac { 13074a54e071SToby Isaac PetscInt p; 13084a54e071SToby Isaac PetscBool valid = PETSC_TRUE; 13094a54e071SToby Isaac for (p = 0; p < newParentSize; ++p) { 13104a54e071SToby Isaac if (newParents[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);} 13114a54e071SToby Isaac } 13124a54e071SToby Isaac if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map"); 13134a54e071SToby Isaac } 13144a54e071SToby Isaac #endif 1315c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr); 1316a6f36705SMichael Lange if (flg) { 1317a6f36705SMichael Lange ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr); 1318a6f36705SMichael Lange ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1319a6f36705SMichael Lange ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr); 1320a6f36705SMichael Lange ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1321a6f36705SMichael Lange ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr); 1322a6f36705SMichael Lange } 1323a6f36705SMichael Lange ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr); 1324a6f36705SMichael Lange ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr); 1325a6f36705SMichael Lange ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr); 1326a6f36705SMichael Lange ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr); 1327a6f36705SMichael Lange } 132815078cd4SMichael Lange pmesh->useAnchors = mesh->useAnchors; 1329a6f36705SMichael Lange PetscFunctionReturn(0); 1330a6f36705SMichael Lange } 13310df0e737SMichael Lange 133224cc2ca5SMatthew G. Knepley PETSC_UNUSED static PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel) 13334eca1733SMichael Lange { 13344eca1733SMichael Lange DM_Plex *mesh = (DM_Plex*) dm->data; 13354eca1733SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 13369852e123SBarry Smith PetscMPIInt rank, size; 13374eca1733SMichael Lange MPI_Comm comm; 13384eca1733SMichael Lange PetscErrorCode ierr; 13394eca1733SMichael Lange 13404eca1733SMichael Lange PetscFunctionBegin; 13414eca1733SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13424eca1733SMichael Lange PetscValidPointer(dmParallel,7); 13434eca1733SMichael Lange 13444eca1733SMichael Lange /* Create point SF for parallel mesh */ 13454eca1733SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 13464eca1733SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 13474eca1733SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 13489852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 13494eca1733SMichael Lange { 13504eca1733SMichael Lange const PetscInt *leaves; 13514eca1733SMichael Lange PetscSFNode *remotePoints, *rowners, *lowners; 13524eca1733SMichael Lange PetscInt numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints; 13534eca1733SMichael Lange PetscInt pStart, pEnd; 13544eca1733SMichael Lange 13554eca1733SMichael Lange ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 13564eca1733SMichael Lange ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr); 13574eca1733SMichael Lange ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr); 13584eca1733SMichael Lange for (p=0; p<numRoots; p++) { 13594eca1733SMichael Lange rowners[p].rank = -1; 13604eca1733SMichael Lange rowners[p].index = -1; 13614eca1733SMichael Lange } 13624eca1733SMichael Lange ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 13634eca1733SMichael Lange ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 13644eca1733SMichael Lange for (p = 0; p < numLeaves; ++p) { 13654eca1733SMichael Lange if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */ 13664eca1733SMichael Lange lowners[p].rank = rank; 13674eca1733SMichael Lange lowners[p].index = leaves ? leaves[p] : p; 13684eca1733SMichael Lange } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */ 13694eca1733SMichael Lange lowners[p].rank = -2; 13704eca1733SMichael Lange lowners[p].index = -2; 13714eca1733SMichael Lange } 13724eca1733SMichael Lange } 13734eca1733SMichael Lange for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */ 13744eca1733SMichael Lange rowners[p].rank = -3; 13754eca1733SMichael Lange rowners[p].index = -3; 13764eca1733SMichael Lange } 13774eca1733SMichael Lange ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 13784eca1733SMichael Lange ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 13794eca1733SMichael Lange ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 13804eca1733SMichael Lange ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 13814eca1733SMichael Lange for (p = 0; p < numLeaves; ++p) { 13824eca1733SMichael Lange if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed"); 13834eca1733SMichael Lange if (lowners[p].rank != rank) ++numGhostPoints; 13844eca1733SMichael Lange } 13854eca1733SMichael Lange ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr); 13864eca1733SMichael Lange ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr); 13874eca1733SMichael Lange for (p = 0, gp = 0; p < numLeaves; ++p) { 13884eca1733SMichael Lange if (lowners[p].rank != rank) { 13894eca1733SMichael Lange ghostPoints[gp] = leaves ? leaves[p] : p; 13904eca1733SMichael Lange remotePoints[gp].rank = lowners[p].rank; 13914eca1733SMichael Lange remotePoints[gp].index = lowners[p].index; 13924eca1733SMichael Lange ++gp; 13934eca1733SMichael Lange } 13944eca1733SMichael Lange } 13954eca1733SMichael Lange ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr); 13964eca1733SMichael Lange ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 13974eca1733SMichael Lange ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr); 13984eca1733SMichael Lange } 13994eca1733SMichael Lange pmesh->useCone = mesh->useCone; 14004eca1733SMichael Lange pmesh->useClosure = mesh->useClosure; 14014eca1733SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 14024eca1733SMichael Lange PetscFunctionReturn(0); 14034eca1733SMichael Lange } 14044eca1733SMichael Lange 1405f5bf2dbfSMichael Lange /*@C 1406f5bf2dbfSMichael Lange DMPlexDerivePointSF - Build a point SF from an SF describing a point migration 1407f5bf2dbfSMichael Lange 1408f5bf2dbfSMichael Lange Input Parameter: 1409f5bf2dbfSMichael Lange + dm - The source DMPlex object 1410f5bf2dbfSMichael Lange . migrationSF - The star forest that describes the parallel point remapping 14111627f6ccSMichael Lange . ownership - Flag causing a vote to determine point ownership 1412f5bf2dbfSMichael Lange 1413f5bf2dbfSMichael Lange Output Parameter: 1414f5bf2dbfSMichael Lange - pointSF - The star forest describing the point overlap in the remapped DM 1415f5bf2dbfSMichael Lange 1416f5bf2dbfSMichael Lange Level: developer 1417f5bf2dbfSMichael Lange 1418f5bf2dbfSMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap() 1419f5bf2dbfSMichael Lange @*/ 1420f5bf2dbfSMichael Lange PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF) 1421f5bf2dbfSMichael Lange { 1422f5bf2dbfSMichael Lange PetscMPIInt rank; 14231627f6ccSMichael Lange PetscInt p, nroots, nleaves, idx, npointLeaves; 1424f5bf2dbfSMichael Lange PetscInt *pointLocal; 1425f5bf2dbfSMichael Lange const PetscInt *leaves; 1426f5bf2dbfSMichael Lange const PetscSFNode *roots; 1427f5bf2dbfSMichael Lange PetscSFNode *rootNodes, *leafNodes, *pointRemote; 1428f5bf2dbfSMichael Lange PetscErrorCode ierr; 1429f5bf2dbfSMichael Lange 1430f5bf2dbfSMichael Lange PetscFunctionBegin; 1431f5bf2dbfSMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1432f5bf2dbfSMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1433f5bf2dbfSMichael Lange 1434f5bf2dbfSMichael Lange ierr = PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots);CHKERRQ(ierr); 1435f5bf2dbfSMichael Lange ierr = PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes);CHKERRQ(ierr); 1436f5bf2dbfSMichael Lange if (ownership) { 14371627f6ccSMichael Lange /* Point ownership vote: Process with highest rank ownes shared points */ 1438f5bf2dbfSMichael Lange for (p = 0; p < nleaves; ++p) { 1439f5bf2dbfSMichael Lange /* Either put in a bid or we know we own it */ 1440f5bf2dbfSMichael Lange leafNodes[p].rank = rank; 144143f7d02bSMichael Lange leafNodes[p].index = p; 1442f5bf2dbfSMichael Lange } 1443f5bf2dbfSMichael Lange for (p = 0; p < nroots; p++) { 14441627f6ccSMichael Lange /* Root must not participate in the reduction, flag so that MAXLOC does not use */ 1445f5bf2dbfSMichael Lange rootNodes[p].rank = -3; 1446f5bf2dbfSMichael Lange rootNodes[p].index = -3; 1447f5bf2dbfSMichael Lange } 1448f5bf2dbfSMichael Lange ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr); 1449f5bf2dbfSMichael Lange ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr); 1450f5bf2dbfSMichael Lange } else { 1451f5bf2dbfSMichael Lange for (p = 0; p < nroots; p++) { 1452f5bf2dbfSMichael Lange rootNodes[p].index = -1; 1453f5bf2dbfSMichael Lange rootNodes[p].rank = rank; 1454f5bf2dbfSMichael Lange }; 1455f5bf2dbfSMichael Lange for (p = 0; p < nleaves; p++) { 1456f5bf2dbfSMichael Lange /* Write new local id into old location */ 1457f5bf2dbfSMichael Lange if (roots[p].rank == rank) { 14586462276dSToby Isaac rootNodes[roots[p].index].index = leaves ? leaves[p] : p; 1459f5bf2dbfSMichael Lange } 1460f5bf2dbfSMichael Lange } 1461f5bf2dbfSMichael Lange } 1462f5bf2dbfSMichael Lange ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr); 1463f5bf2dbfSMichael Lange ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr); 1464f5bf2dbfSMichael Lange 1465f5bf2dbfSMichael Lange for (npointLeaves = 0, p = 0; p < nleaves; p++) {if (leafNodes[p].rank != rank) npointLeaves++;} 14661627f6ccSMichael Lange ierr = PetscMalloc1(npointLeaves, &pointLocal);CHKERRQ(ierr); 14671627f6ccSMichael Lange ierr = PetscMalloc1(npointLeaves, &pointRemote);CHKERRQ(ierr); 1468f5bf2dbfSMichael Lange for (idx = 0, p = 0; p < nleaves; p++) { 1469f5bf2dbfSMichael Lange if (leafNodes[p].rank != rank) { 1470f5bf2dbfSMichael Lange pointLocal[idx] = p; 1471f5bf2dbfSMichael Lange pointRemote[idx] = leafNodes[p]; 1472f5bf2dbfSMichael Lange idx++; 1473f5bf2dbfSMichael Lange } 1474f5bf2dbfSMichael Lange } 1475f5bf2dbfSMichael Lange ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF);CHKERRQ(ierr); 14761627f6ccSMichael Lange ierr = PetscSFSetFromOptions(*pointSF);CHKERRQ(ierr); 1477f5bf2dbfSMichael Lange ierr = PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER);CHKERRQ(ierr); 1478f5bf2dbfSMichael Lange ierr = PetscFree2(rootNodes, leafNodes);CHKERRQ(ierr); 1479f5bf2dbfSMichael Lange PetscFunctionReturn(0); 1480f5bf2dbfSMichael Lange } 1481f5bf2dbfSMichael Lange 148215078cd4SMichael Lange /*@C 148315078cd4SMichael Lange DMPlexMigrate - Migrates internal DM data over the supplied star forest 148415078cd4SMichael Lange 148515078cd4SMichael Lange Input Parameter: 148615078cd4SMichael Lange + dm - The source DMPlex object 14871627f6ccSMichael Lange . sf - The star forest communication context describing the migration pattern 148815078cd4SMichael Lange 148915078cd4SMichael Lange Output Parameter: 149015078cd4SMichael Lange - targetDM - The target DMPlex object 149115078cd4SMichael Lange 1492b9f40539SMichael Lange Level: intermediate 149315078cd4SMichael Lange 149415078cd4SMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap() 149515078cd4SMichael Lange @*/ 1496b9f40539SMichael Lange PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM) 149715078cd4SMichael Lange { 1498b9f40539SMichael Lange MPI_Comm comm; 1499b9f40539SMichael Lange PetscInt dim, nroots; 1500b9f40539SMichael Lange PetscSF sfPoint; 150115078cd4SMichael Lange ISLocalToGlobalMapping ltogMigration; 1502ac04eaf7SToby Isaac ISLocalToGlobalMapping ltogOriginal = NULL; 1503b9f40539SMichael Lange PetscBool flg; 150415078cd4SMichael Lange PetscErrorCode ierr; 150515078cd4SMichael Lange 150615078cd4SMichael Lange PetscFunctionBegin; 150715078cd4SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15081b858b30SMichael Lange ierr = PetscLogEventBegin(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr); 1509b9f40539SMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 151015078cd4SMichael Lange ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 151115078cd4SMichael Lange ierr = DMSetDimension(targetDM, dim);CHKERRQ(ierr); 151215078cd4SMichael Lange 1513bfb0467fSMichael Lange /* Check for a one-to-all distribution pattern */ 1514b9f40539SMichael Lange ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1515b9f40539SMichael Lange ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1516bfb0467fSMichael Lange if (nroots >= 0) { 1517b9f40539SMichael Lange IS isOriginal; 1518ac04eaf7SToby Isaac PetscInt n, size, nleaves; 1519ac04eaf7SToby Isaac PetscInt *numbering_orig, *numbering_new; 1520b9f40539SMichael Lange /* Get the original point numbering */ 1521b9f40539SMichael Lange ierr = DMPlexCreatePointNumbering(dm, &isOriginal);CHKERRQ(ierr); 1522b9f40539SMichael Lange ierr = ISLocalToGlobalMappingCreateIS(isOriginal, <ogOriginal);CHKERRQ(ierr); 1523b9f40539SMichael Lange ierr = ISLocalToGlobalMappingGetSize(ltogOriginal, &size);CHKERRQ(ierr); 1524b9f40539SMichael Lange ierr = ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr); 1525b9f40539SMichael Lange /* Convert to positive global numbers */ 1526b9f40539SMichael Lange for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);} 1527b9f40539SMichael Lange /* Derive the new local-to-global mapping from the old one */ 1528b9f40539SMichael Lange ierr = PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr); 1529b9f40539SMichael Lange ierr = PetscMalloc1(nleaves, &numbering_new);CHKERRQ(ierr); 1530b9f40539SMichael Lange ierr = PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr); 1531b9f40539SMichael Lange ierr = PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr); 1532b9f40539SMichael Lange ierr = ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, <ogMigration);CHKERRQ(ierr); 1533b9f40539SMichael Lange ierr = ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr); 1534b9f40539SMichael Lange ierr = ISDestroy(&isOriginal);CHKERRQ(ierr); 153515078cd4SMichael Lange } else { 1536bfb0467fSMichael Lange /* One-to-all distribution pattern: We can derive LToG from SF */ 153715078cd4SMichael Lange ierr = ISLocalToGlobalMappingCreateSF(sf, 0, <ogMigration);CHKERRQ(ierr); 153815078cd4SMichael Lange } 1539c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 1540b9f40539SMichael Lange if (flg) { 1541b9f40539SMichael Lange ierr = PetscPrintf(comm, "Point renumbering for DM migration:\n");CHKERRQ(ierr); 1542b9f40539SMichael Lange ierr = ISLocalToGlobalMappingView(ltogMigration, NULL);CHKERRQ(ierr); 1543b9f40539SMichael Lange } 154415078cd4SMichael Lange /* Migrate DM data to target DM */ 1545ac04eaf7SToby Isaac ierr = DMPlexDistributeCones(dm, sf, ltogOriginal, ltogMigration, targetDM);CHKERRQ(ierr); 154615078cd4SMichael Lange ierr = DMPlexDistributeLabels(dm, sf, targetDM);CHKERRQ(ierr); 15470dc1530dSSander Arens ierr = DMPlexDistributeCoordinates(dm, sf, targetDM);CHKERRQ(ierr); 154815078cd4SMichael Lange ierr = DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr); 154908633170SToby Isaac ierr = DMPlexDistributeSetupTree(dm, sf, ltogOriginal, ltogMigration, targetDM);CHKERRQ(ierr); 1550ac04eaf7SToby Isaac ierr = ISLocalToGlobalMappingDestroy(<ogOriginal);CHKERRQ(ierr); 1551302440fdSBarry Smith ierr = ISLocalToGlobalMappingDestroy(<ogMigration);CHKERRQ(ierr); 15521b858b30SMichael Lange ierr = PetscLogEventEnd(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr); 155315078cd4SMichael Lange PetscFunctionReturn(0); 155415078cd4SMichael Lange } 155515078cd4SMichael Lange 15563b8f15a2SMatthew G. Knepley /*@C 155770034214SMatthew G. Knepley DMPlexDistribute - Distributes the mesh and any associated sections. 155870034214SMatthew G. Knepley 155970034214SMatthew G. Knepley Not Collective 156070034214SMatthew G. Knepley 156170034214SMatthew G. Knepley Input Parameter: 156270034214SMatthew G. Knepley + dm - The original DMPlex object 156370034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default 156470034214SMatthew G. Knepley 156570034214SMatthew G. Knepley Output Parameter: 156670034214SMatthew G. Knepley + sf - The PetscSF used for point distribution 156770034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL 156870034214SMatthew G. Knepley 156970034214SMatthew G. Knepley Note: If the mesh was not distributed, the return value is NULL. 157070034214SMatthew G. Knepley 157193a1dc33SMatthew G. Knepley The user can control the definition of adjacency for the mesh using DMPlexSetAdjacencyUseCone() and 157270034214SMatthew G. Knepley DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 157370034214SMatthew G. Knepley representation on the mesh. 157470034214SMatthew G. Knepley 157570034214SMatthew G. Knepley Level: intermediate 157670034214SMatthew G. Knepley 157770034214SMatthew G. Knepley .keywords: mesh, elements 157870034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 157970034214SMatthew G. Knepley @*/ 158080cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel) 158170034214SMatthew G. Knepley { 158270034214SMatthew G. Knepley MPI_Comm comm; 158315078cd4SMichael Lange PetscPartitioner partitioner; 1584f8987ae8SMichael Lange IS cellPart; 1585f8987ae8SMichael Lange PetscSection cellPartSection; 1586cf86098cSMatthew G. Knepley DM dmCoord; 1587f8987ae8SMichael Lange DMLabel lblPartition, lblMigration; 158843f7d02bSMichael Lange PetscSF sfProcess, sfMigration, sfStratified, sfPoint; 158970034214SMatthew G. Knepley PetscBool flg; 15909852e123SBarry Smith PetscMPIInt rank, size, p; 159170034214SMatthew G. Knepley PetscErrorCode ierr; 159270034214SMatthew G. Knepley 159370034214SMatthew G. Knepley PetscFunctionBegin; 159470034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 159570034214SMatthew G. Knepley if (sf) PetscValidPointer(sf,4); 159670034214SMatthew G. Knepley PetscValidPointer(dmParallel,5); 159770034214SMatthew G. Knepley 159870034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 159970034214SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 160070034214SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 16019852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 160270034214SMatthew G. Knepley 1603*a6ba33f0SVaclav Hapla if (sf) *sf = NULL; 160470034214SMatthew G. Knepley *dmParallel = NULL; 16059852e123SBarry Smith if (size == 1) { 16065f3267c8SKoos Huijssen ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 16075f3267c8SKoos Huijssen PetscFunctionReturn(0); 16085f3267c8SKoos Huijssen } 160970034214SMatthew G. Knepley 161015078cd4SMichael Lange /* Create cell partition */ 161177623264SMatthew G. Knepley ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 161277623264SMatthew G. Knepley ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr); 161315078cd4SMichael Lange ierr = DMPlexGetPartitioner(dm, &partitioner);CHKERRQ(ierr); 161415078cd4SMichael Lange ierr = PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);CHKERRQ(ierr); 1615f8987ae8SMichael Lange { 1616f8987ae8SMichael Lange /* Convert partition to DMLabel */ 1617f8987ae8SMichael Lange PetscInt proc, pStart, pEnd, npoints, poffset; 1618f8987ae8SMichael Lange const PetscInt *points; 1619f8987ae8SMichael Lange ierr = DMLabelCreate("Point Partition", &lblPartition);CHKERRQ(ierr); 1620f8987ae8SMichael Lange ierr = ISGetIndices(cellPart, &points);CHKERRQ(ierr); 1621f8987ae8SMichael Lange ierr = PetscSectionGetChart(cellPartSection, &pStart, &pEnd);CHKERRQ(ierr); 1622f8987ae8SMichael Lange for (proc = pStart; proc < pEnd; proc++) { 1623f8987ae8SMichael Lange ierr = PetscSectionGetDof(cellPartSection, proc, &npoints);CHKERRQ(ierr); 1624f8987ae8SMichael Lange ierr = PetscSectionGetOffset(cellPartSection, proc, &poffset);CHKERRQ(ierr); 1625f8987ae8SMichael Lange for (p = poffset; p < poffset+npoints; p++) { 1626f8987ae8SMichael Lange ierr = DMLabelSetValue(lblPartition, points[p], proc);CHKERRQ(ierr); 162770034214SMatthew G. Knepley } 1628f8987ae8SMichael Lange } 1629f8987ae8SMichael Lange ierr = ISRestoreIndices(cellPart, &points);CHKERRQ(ierr); 1630f8987ae8SMichael Lange } 1631f8987ae8SMichael Lange ierr = DMPlexPartitionLabelClosure(dm, lblPartition);CHKERRQ(ierr); 1632f8987ae8SMichael Lange { 1633f8987ae8SMichael Lange /* Build a global process SF */ 1634f8987ae8SMichael Lange PetscSFNode *remoteProc; 16359852e123SBarry Smith ierr = PetscMalloc1(size, &remoteProc);CHKERRQ(ierr); 16369852e123SBarry Smith for (p = 0; p < size; ++p) { 1637f8987ae8SMichael Lange remoteProc[p].rank = p; 1638f8987ae8SMichael Lange remoteProc[p].index = rank; 1639f8987ae8SMichael Lange } 1640f8987ae8SMichael Lange ierr = PetscSFCreate(comm, &sfProcess);CHKERRQ(ierr); 1641f8987ae8SMichael Lange ierr = PetscObjectSetName((PetscObject) sfProcess, "Process SF");CHKERRQ(ierr); 16429852e123SBarry Smith ierr = PetscSFSetGraph(sfProcess, size, size, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr); 1643f8987ae8SMichael Lange } 1644f8987ae8SMichael Lange ierr = DMLabelCreate("Point migration", &lblMigration);CHKERRQ(ierr); 1645f8987ae8SMichael Lange ierr = DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);CHKERRQ(ierr); 1646302440fdSBarry Smith ierr = DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration);CHKERRQ(ierr); 164743f7d02bSMichael Lange /* Stratify the SF in case we are migrating an already parallel plex */ 164843f7d02bSMichael Lange ierr = DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified);CHKERRQ(ierr); 164943f7d02bSMichael Lange ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr); 165043f7d02bSMichael Lange sfMigration = sfStratified; 1651f8987ae8SMichael Lange ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1652c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 165370034214SMatthew G. Knepley if (flg) { 1654f8987ae8SMichael Lange ierr = DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1655f8987ae8SMichael Lange ierr = PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 165670034214SMatthew G. Knepley } 1657f8987ae8SMichael Lange 165815078cd4SMichael Lange /* Create non-overlapping parallel DM and migrate internal data */ 165970034214SMatthew G. Knepley ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr); 166070034214SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr); 1661b9f40539SMichael Lange ierr = DMPlexMigrate(dm, sfMigration, *dmParallel);CHKERRQ(ierr); 166270034214SMatthew G. Knepley 1663a157612eSMichael Lange /* Build the point SF without overlap */ 16641627f6ccSMichael Lange ierr = DMPlexCreatePointSF(*dmParallel, sfMigration, PETSC_TRUE, &sfPoint);CHKERRQ(ierr); 1665f5bf2dbfSMichael Lange ierr = DMSetPointSF(*dmParallel, sfPoint);CHKERRQ(ierr); 1666cf86098cSMatthew G. Knepley ierr = DMGetCoordinateDM(*dmParallel, &dmCoord);CHKERRQ(ierr); 1667cf86098cSMatthew G. Knepley if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);} 1668f5bf2dbfSMichael Lange if (flg) {ierr = PetscSFView(sfPoint, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 166970034214SMatthew G. Knepley 1670a157612eSMichael Lange if (overlap > 0) { 167115078cd4SMichael Lange DM dmOverlap; 167215078cd4SMichael Lange PetscInt nroots, nleaves; 167315078cd4SMichael Lange PetscSFNode *newRemote; 167415078cd4SMichael Lange const PetscSFNode *oldRemote; 167515078cd4SMichael Lange PetscSF sfOverlap, sfOverlapPoint; 1676a157612eSMichael Lange /* Add the partition overlap to the distributed DM */ 1677b9f40539SMichael Lange ierr = DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);CHKERRQ(ierr); 1678a157612eSMichael Lange ierr = DMDestroy(dmParallel);CHKERRQ(ierr); 1679a157612eSMichael Lange *dmParallel = dmOverlap; 1680c389ea9fSToby Isaac if (flg) { 16813d822a50SMichael Lange ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr); 168215078cd4SMichael Lange ierr = PetscSFView(sfOverlap, NULL);CHKERRQ(ierr); 1683c389ea9fSToby Isaac } 168443331d4aSMichael Lange 1685f8987ae8SMichael Lange /* Re-map the migration SF to establish the full migration pattern */ 1686f8987ae8SMichael Lange ierr = PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr); 168715078cd4SMichael Lange ierr = PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr); 168843331d4aSMichael Lange ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr); 168915078cd4SMichael Lange ierr = PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr); 169015078cd4SMichael Lange ierr = PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr); 169115078cd4SMichael Lange ierr = PetscSFCreate(comm, &sfOverlapPoint);CHKERRQ(ierr); 169215078cd4SMichael Lange ierr = PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr); 169315078cd4SMichael Lange ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr); 1694f8987ae8SMichael Lange ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr); 169515078cd4SMichael Lange sfMigration = sfOverlapPoint; 1696c389ea9fSToby Isaac } 1697bf5244c7SMatthew G. Knepley /* Cleanup Partition */ 1698f8987ae8SMichael Lange ierr = PetscSFDestroy(&sfProcess);CHKERRQ(ierr); 1699f8987ae8SMichael Lange ierr = DMLabelDestroy(&lblPartition);CHKERRQ(ierr); 1700f8987ae8SMichael Lange ierr = DMLabelDestroy(&lblMigration);CHKERRQ(ierr); 17014eca1733SMichael Lange ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); 17024eca1733SMichael Lange ierr = ISDestroy(&cellPart);CHKERRQ(ierr); 1703721cbd76SMatthew G. Knepley /* Copy BC */ 1704a6ba4734SToby Isaac ierr = DMCopyBoundary(dm, *dmParallel);CHKERRQ(ierr); 170566fe0bfeSMatthew G. Knepley /* Create sfNatural */ 170666fe0bfeSMatthew G. Knepley if (dm->useNatural) { 170766fe0bfeSMatthew G. Knepley PetscSection section; 170866fe0bfeSMatthew G. Knepley 170966fe0bfeSMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 17100c1f158bSMatthew G. Knepley ierr = DMPlexCreateGlobalToNaturalSF(*dmParallel, section, sfMigration, &(*dmParallel)->sfNatural);CHKERRQ(ierr); 171166fe0bfeSMatthew G. Knepley } 1712721cbd76SMatthew G. Knepley /* Cleanup */ 1713f8987ae8SMichael Lange if (sf) {*sf = sfMigration;} 1714f8987ae8SMichael Lange else {ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);} 1715f5bf2dbfSMichael Lange ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr); 171670034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 171770034214SMatthew G. Knepley PetscFunctionReturn(0); 171870034214SMatthew G. Knepley } 1719a157612eSMichael Lange 1720a157612eSMichael Lange /*@C 172124cc2ca5SMatthew G. Knepley DMPlexDistributeOverlap - Add partition overlap to a distributed non-overlapping DM. 1722a157612eSMichael Lange 1723a157612eSMichael Lange Not Collective 1724a157612eSMichael Lange 1725a157612eSMichael Lange Input Parameter: 1726a157612eSMichael Lange + dm - The non-overlapping distrbuted DMPlex object 1727a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default 1728a157612eSMichael Lange 1729a157612eSMichael Lange Output Parameter: 1730a157612eSMichael Lange + sf - The PetscSF used for point distribution 1731a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL 1732a157612eSMichael Lange 1733a157612eSMichael Lange Note: If the mesh was not distributed, the return value is NULL. 1734a157612eSMichael Lange 1735a157612eSMichael Lange The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 1736a157612eSMichael Lange DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 1737a157612eSMichael Lange representation on the mesh. 1738a157612eSMichael Lange 1739a157612eSMichael Lange Level: intermediate 1740a157612eSMichael Lange 1741a157612eSMichael Lange .keywords: mesh, elements 1742a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 1743a157612eSMichael Lange @*/ 1744b9f40539SMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap) 1745a157612eSMichael Lange { 1746a157612eSMichael Lange MPI_Comm comm; 17473567eaeeSMatthew G. Knepley PetscMPIInt size, rank; 1748a157612eSMichael Lange PetscSection rootSection, leafSection; 1749a157612eSMichael Lange IS rootrank, leafrank; 1750cf86098cSMatthew G. Knepley DM dmCoord; 1751a9f1d5b2SMichael Lange DMLabel lblOverlap; 1752f5bf2dbfSMichael Lange PetscSF sfOverlap, sfStratified, sfPoint; 1753a157612eSMichael Lange PetscErrorCode ierr; 1754a157612eSMichael Lange 1755a157612eSMichael Lange PetscFunctionBegin; 1756a157612eSMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1757a157612eSMichael Lange if (sf) PetscValidPointer(sf, 3); 1758a157612eSMichael Lange PetscValidPointer(dmOverlap, 4); 1759a157612eSMichael Lange 1760a157612eSMichael Lange ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 17613567eaeeSMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1762a157612eSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 17633567eaeeSMatthew G. Knepley if (size == 1) {*dmOverlap = NULL; PetscFunctionReturn(0);} 17643567eaeeSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr); 1765a157612eSMichael Lange 1766a157612eSMichael Lange /* Compute point overlap with neighbouring processes on the distributed DM */ 1767a157612eSMichael Lange ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1768a157612eSMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 1769a157612eSMichael Lange ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); 1770a157612eSMichael Lange ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr); 1771a9f1d5b2SMichael Lange ierr = DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap);CHKERRQ(ierr); 1772a9f1d5b2SMichael Lange /* Convert overlap label to stratified migration SF */ 1773a9f1d5b2SMichael Lange ierr = DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap);CHKERRQ(ierr); 1774a9f1d5b2SMichael Lange ierr = DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified);CHKERRQ(ierr); 1775a9f1d5b2SMichael Lange ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr); 1776a9f1d5b2SMichael Lange sfOverlap = sfStratified; 1777a9f1d5b2SMichael Lange ierr = PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF");CHKERRQ(ierr); 1778a9f1d5b2SMichael Lange ierr = PetscSFSetFromOptions(sfOverlap);CHKERRQ(ierr); 1779a9f1d5b2SMichael Lange 178015fff7beSMatthew G. Knepley ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 178115fff7beSMatthew G. Knepley ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 178215fff7beSMatthew G. Knepley ierr = ISDestroy(&rootrank);CHKERRQ(ierr); 178315fff7beSMatthew G. Knepley ierr = ISDestroy(&leafrank);CHKERRQ(ierr); 1784a157612eSMichael Lange ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1785a157612eSMichael Lange 1786a157612eSMichael Lange /* Build the overlapping DM */ 1787a157612eSMichael Lange ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr); 1788a157612eSMichael Lange ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr); 1789a9f1d5b2SMichael Lange ierr = DMPlexMigrate(dm, sfOverlap, *dmOverlap);CHKERRQ(ierr); 1790f5bf2dbfSMichael Lange /* Build the new point SF */ 17911627f6ccSMichael Lange ierr = DMPlexCreatePointSF(*dmOverlap, sfOverlap, PETSC_FALSE, &sfPoint);CHKERRQ(ierr); 1792f5bf2dbfSMichael Lange ierr = DMSetPointSF(*dmOverlap, sfPoint);CHKERRQ(ierr); 1793cf86098cSMatthew G. Knepley ierr = DMGetCoordinateDM(*dmOverlap, &dmCoord);CHKERRQ(ierr); 1794cf86098cSMatthew G. Knepley if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);} 1795f5bf2dbfSMichael Lange ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr); 1796a157612eSMichael Lange /* Cleanup overlap partition */ 1797a9f1d5b2SMichael Lange ierr = DMLabelDestroy(&lblOverlap);CHKERRQ(ierr); 1798a9f1d5b2SMichael Lange if (sf) *sf = sfOverlap; 1799a9f1d5b2SMichael Lange else {ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);} 18001b858b30SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr); 1801a157612eSMichael Lange PetscFunctionReturn(0); 1802a157612eSMichael Lange } 18036462276dSToby Isaac 18046462276dSToby Isaac /*@C 18056462276dSToby Isaac DMPlexGetGatherDM - Get a copy of the DMPlex that gathers all points on the 18066462276dSToby Isaac root process of the original's communicator. 18076462276dSToby Isaac 18086462276dSToby Isaac Input Parameters: 18096462276dSToby Isaac . dm - the original DMPlex object 18106462276dSToby Isaac 18116462276dSToby Isaac Output Parameters: 18126462276dSToby Isaac . gatherMesh - the gathered DM object, or NULL 18136462276dSToby Isaac 18146462276dSToby Isaac Level: intermediate 18156462276dSToby Isaac 18166462276dSToby Isaac .keywords: mesh 18176462276dSToby Isaac .seealso: DMPlexDistribute(), DMPlexGetRedundantDM() 18186462276dSToby Isaac @*/ 18196462276dSToby Isaac PetscErrorCode DMPlexGetGatherDM(DM dm, DM * gatherMesh) 18206462276dSToby Isaac { 18216462276dSToby Isaac MPI_Comm comm; 18226462276dSToby Isaac PetscMPIInt size; 18236462276dSToby Isaac PetscPartitioner oldPart, gatherPart; 18246462276dSToby Isaac PetscErrorCode ierr; 18256462276dSToby Isaac 18266462276dSToby Isaac PetscFunctionBegin; 18276462276dSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 18286462276dSToby Isaac comm = PetscObjectComm((PetscObject)dm); 18296462276dSToby Isaac ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 18306462276dSToby Isaac *gatherMesh = NULL; 18316462276dSToby Isaac if (size == 1) PetscFunctionReturn(0); 18326462276dSToby Isaac ierr = DMPlexGetPartitioner(dm,&oldPart);CHKERRQ(ierr); 18336462276dSToby Isaac ierr = PetscObjectReference((PetscObject)oldPart);CHKERRQ(ierr); 18346462276dSToby Isaac ierr = PetscPartitionerCreate(comm,&gatherPart);CHKERRQ(ierr); 18356462276dSToby Isaac ierr = PetscPartitionerSetType(gatherPart,PETSCPARTITIONERGATHER);CHKERRQ(ierr); 18366462276dSToby Isaac ierr = DMPlexSetPartitioner(dm,gatherPart);CHKERRQ(ierr); 18376462276dSToby Isaac ierr = DMPlexDistribute(dm,0,NULL,gatherMesh);CHKERRQ(ierr); 18386462276dSToby Isaac ierr = DMPlexSetPartitioner(dm,oldPart);CHKERRQ(ierr); 18396462276dSToby Isaac ierr = PetscPartitionerDestroy(&gatherPart);CHKERRQ(ierr); 18406462276dSToby Isaac ierr = PetscPartitionerDestroy(&oldPart);CHKERRQ(ierr); 18416462276dSToby Isaac PetscFunctionReturn(0); 18426462276dSToby Isaac } 18436462276dSToby Isaac 18446462276dSToby Isaac /*@C 18456462276dSToby Isaac DMPlexGetRedundantDM - Get a copy of the DMPlex that is completely copied on each process. 18466462276dSToby Isaac 18476462276dSToby Isaac Input Parameters: 18486462276dSToby Isaac . dm - the original DMPlex object 18496462276dSToby Isaac 18506462276dSToby Isaac Output Parameters: 18516462276dSToby Isaac . redundantMesh - the redundant DM object, or NULL 18526462276dSToby Isaac 18536462276dSToby Isaac Level: intermediate 18546462276dSToby Isaac 18556462276dSToby Isaac .keywords: mesh 18566462276dSToby Isaac .seealso: DMPlexDistribute(), DMPlexGetGatherDM() 18576462276dSToby Isaac @*/ 18586462276dSToby Isaac PetscErrorCode DMPlexGetRedundantDM(DM dm, DM * redundantMesh) 18596462276dSToby Isaac { 18606462276dSToby Isaac MPI_Comm comm; 18616462276dSToby Isaac PetscMPIInt size, rank; 18626462276dSToby Isaac PetscInt pStart, pEnd, p; 18636462276dSToby Isaac PetscInt numPoints = -1; 18646462276dSToby Isaac PetscSF migrationSF, sfPoint; 18656462276dSToby Isaac DM gatherDM, dmCoord; 18666462276dSToby Isaac PetscSFNode *points; 18676462276dSToby Isaac PetscErrorCode ierr; 18686462276dSToby Isaac 18696462276dSToby Isaac PetscFunctionBegin; 18706462276dSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 18716462276dSToby Isaac comm = PetscObjectComm((PetscObject)dm); 18726462276dSToby Isaac ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 18736462276dSToby Isaac *redundantMesh = NULL; 187468dbc166SMatthew G. Knepley if (size == 1) { 187568dbc166SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 187668dbc166SMatthew G. Knepley *redundantMesh = dm; 187768dbc166SMatthew G. Knepley PetscFunctionReturn(0); 187868dbc166SMatthew G. Knepley } 18796462276dSToby Isaac ierr = DMPlexGetGatherDM(dm,&gatherDM);CHKERRQ(ierr); 18806462276dSToby Isaac if (!gatherDM) PetscFunctionReturn(0); 18816462276dSToby Isaac ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 18826462276dSToby Isaac ierr = DMPlexGetChart(gatherDM,&pStart,&pEnd);CHKERRQ(ierr); 18836462276dSToby Isaac numPoints = pEnd - pStart; 18846462276dSToby Isaac ierr = MPI_Bcast(&numPoints,1,MPIU_INT,0,comm);CHKERRQ(ierr); 18856462276dSToby Isaac ierr = PetscMalloc1(numPoints,&points);CHKERRQ(ierr); 18866462276dSToby Isaac ierr = PetscSFCreate(comm,&migrationSF);CHKERRQ(ierr); 18876462276dSToby Isaac for (p = 0; p < numPoints; p++) { 18886462276dSToby Isaac points[p].index = p; 18896462276dSToby Isaac points[p].rank = 0; 18906462276dSToby Isaac } 18916462276dSToby Isaac ierr = PetscSFSetGraph(migrationSF,pEnd-pStart,numPoints,NULL,PETSC_OWN_POINTER,points,PETSC_OWN_POINTER);CHKERRQ(ierr); 18926462276dSToby Isaac ierr = DMPlexCreate(comm, redundantMesh);CHKERRQ(ierr); 18936462276dSToby Isaac ierr = PetscObjectSetName((PetscObject) *redundantMesh, "Redundant Mesh");CHKERRQ(ierr); 18946462276dSToby Isaac ierr = DMPlexMigrate(gatherDM, migrationSF, *redundantMesh);CHKERRQ(ierr); 18956462276dSToby Isaac ierr = DMPlexCreatePointSF(*redundantMesh, migrationSF, PETSC_FALSE, &sfPoint);CHKERRQ(ierr); 18966462276dSToby Isaac ierr = DMSetPointSF(*redundantMesh, sfPoint);CHKERRQ(ierr); 18976462276dSToby Isaac ierr = DMGetCoordinateDM(*redundantMesh, &dmCoord);CHKERRQ(ierr); 18986462276dSToby Isaac if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);} 18996462276dSToby Isaac ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr); 19006462276dSToby Isaac ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr); 19016462276dSToby Isaac ierr = DMDestroy(&gatherDM);CHKERRQ(ierr); 19026462276dSToby Isaac PetscFunctionReturn(0); 19036462276dSToby Isaac } 1904