170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 270034214SMatthew G. Knepley 370034214SMatthew G. Knepley #undef __FUNCT__ 470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone" 570034214SMatthew G. Knepley /*@ 670034214SMatthew G. Knepley DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first 770034214SMatthew G. Knepley 870034214SMatthew G. Knepley Input Parameters: 970034214SMatthew G. Knepley + dm - The DM object 1070034214SMatthew G. Knepley - useCone - Flag to use the cone first 1170034214SMatthew G. Knepley 1270034214SMatthew G. Knepley Level: intermediate 1370034214SMatthew G. Knepley 1470034214SMatthew G. Knepley Notes: 1570034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 164b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 1770034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 1870034214SMatthew G. Knepley 1970034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 2070034214SMatthew G. Knepley @*/ 2170034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone) 2270034214SMatthew G. Knepley { 2370034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 2470034214SMatthew G. Knepley 2570034214SMatthew G. Knepley PetscFunctionBegin; 2670034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2770034214SMatthew G. Knepley mesh->useCone = useCone; 2870034214SMatthew G. Knepley PetscFunctionReturn(0); 2970034214SMatthew G. Knepley } 3070034214SMatthew G. Knepley 3170034214SMatthew G. Knepley #undef __FUNCT__ 3270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone" 3370034214SMatthew G. Knepley /*@ 3470034214SMatthew G. Knepley DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first 3570034214SMatthew G. Knepley 3670034214SMatthew G. Knepley Input Parameter: 3770034214SMatthew G. Knepley . dm - The DM object 3870034214SMatthew G. Knepley 3970034214SMatthew G. Knepley Output Parameter: 4070034214SMatthew G. Knepley . useCone - Flag to use the cone first 4170034214SMatthew G. Knepley 4270034214SMatthew G. Knepley Level: intermediate 4370034214SMatthew G. Knepley 4470034214SMatthew G. Knepley Notes: 4570034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 464b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 4770034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 4870034214SMatthew G. Knepley 4970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 5070034214SMatthew G. Knepley @*/ 5170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone) 5270034214SMatthew G. Knepley { 5370034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 5470034214SMatthew G. Knepley 5570034214SMatthew G. Knepley PetscFunctionBegin; 5670034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5770034214SMatthew G. Knepley PetscValidIntPointer(useCone, 2); 5870034214SMatthew G. Knepley *useCone = mesh->useCone; 5970034214SMatthew G. Knepley PetscFunctionReturn(0); 6070034214SMatthew G. Knepley } 6170034214SMatthew G. Knepley 6270034214SMatthew G. Knepley #undef __FUNCT__ 6370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure" 6470034214SMatthew G. Knepley /*@ 6570034214SMatthew G. Knepley DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure 6670034214SMatthew G. Knepley 6770034214SMatthew G. Knepley Input Parameters: 6870034214SMatthew G. Knepley + dm - The DM object 6970034214SMatthew G. Knepley - useClosure - Flag to use the closure 7070034214SMatthew G. Knepley 7170034214SMatthew G. Knepley Level: intermediate 7270034214SMatthew G. Knepley 7370034214SMatthew G. Knepley Notes: 7470034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 754b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 7670034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 7770034214SMatthew G. Knepley 7870034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 7970034214SMatthew G. Knepley @*/ 8070034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure) 8170034214SMatthew G. Knepley { 8270034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 8370034214SMatthew G. Knepley 8470034214SMatthew G. Knepley PetscFunctionBegin; 8570034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8670034214SMatthew G. Knepley mesh->useClosure = useClosure; 8770034214SMatthew G. Knepley PetscFunctionReturn(0); 8870034214SMatthew G. Knepley } 8970034214SMatthew G. Knepley 9070034214SMatthew G. Knepley #undef __FUNCT__ 9170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure" 9270034214SMatthew G. Knepley /*@ 9370034214SMatthew G. Knepley DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure 9470034214SMatthew G. Knepley 9570034214SMatthew G. Knepley Input Parameter: 9670034214SMatthew G. Knepley . dm - The DM object 9770034214SMatthew G. Knepley 9870034214SMatthew G. Knepley Output Parameter: 9970034214SMatthew G. Knepley . useClosure - Flag to use the closure 10070034214SMatthew G. Knepley 10170034214SMatthew G. Knepley Level: intermediate 10270034214SMatthew G. Knepley 10370034214SMatthew G. Knepley Notes: 10470034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 1054b6b44bdSMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 10670034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 10770034214SMatthew G. Knepley 10870034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 10970034214SMatthew G. Knepley @*/ 11070034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure) 11170034214SMatthew G. Knepley { 11270034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 11370034214SMatthew G. Knepley 11470034214SMatthew G. Knepley PetscFunctionBegin; 11570034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11670034214SMatthew G. Knepley PetscValidIntPointer(useClosure, 2); 11770034214SMatthew G. Knepley *useClosure = mesh->useClosure; 11870034214SMatthew G. Knepley PetscFunctionReturn(0); 11970034214SMatthew G. Knepley } 12070034214SMatthew G. Knepley 12170034214SMatthew G. Knepley #undef __FUNCT__ 122a17985deSToby Isaac #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors" 1238b0b4c70SToby Isaac /*@ 124a17985deSToby Isaac DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints. 1258b0b4c70SToby Isaac 1268b0b4c70SToby Isaac Input Parameters: 1278b0b4c70SToby Isaac + dm - The DM object 1285b317d89SToby 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. 1298b0b4c70SToby Isaac 1308b0b4c70SToby Isaac Level: intermediate 1318b0b4c70SToby Isaac 132a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors() 1338b0b4c70SToby Isaac @*/ 1345b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors) 1358b0b4c70SToby Isaac { 1368b0b4c70SToby Isaac DM_Plex *mesh = (DM_Plex *) dm->data; 1378b0b4c70SToby Isaac 1388b0b4c70SToby Isaac PetscFunctionBegin; 1398b0b4c70SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1405b317d89SToby Isaac mesh->useAnchors = useAnchors; 1418b0b4c70SToby Isaac PetscFunctionReturn(0); 1428b0b4c70SToby Isaac } 1438b0b4c70SToby Isaac 1448b0b4c70SToby Isaac #undef __FUNCT__ 145a17985deSToby Isaac #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors" 1468b0b4c70SToby Isaac /*@ 147a17985deSToby Isaac DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints. 1488b0b4c70SToby Isaac 1498b0b4c70SToby Isaac Input Parameter: 1508b0b4c70SToby Isaac . dm - The DM object 1518b0b4c70SToby Isaac 1528b0b4c70SToby Isaac Output Parameter: 1535b317d89SToby 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. 1548b0b4c70SToby Isaac 1558b0b4c70SToby Isaac Level: intermediate 1568b0b4c70SToby Isaac 157a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors() 1588b0b4c70SToby Isaac @*/ 1595b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors) 1608b0b4c70SToby Isaac { 1618b0b4c70SToby Isaac DM_Plex *mesh = (DM_Plex *) dm->data; 1628b0b4c70SToby Isaac 1638b0b4c70SToby Isaac PetscFunctionBegin; 1648b0b4c70SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1655b317d89SToby Isaac PetscValidIntPointer(useAnchors, 2); 1665b317d89SToby Isaac *useAnchors = mesh->useAnchors; 1678b0b4c70SToby Isaac PetscFunctionReturn(0); 1688b0b4c70SToby Isaac } 1698b0b4c70SToby Isaac 1708b0b4c70SToby Isaac #undef __FUNCT__ 17170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal" 17270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 17370034214SMatthew G. Knepley { 17470034214SMatthew G. Knepley const PetscInt *cone = NULL; 17570034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, coneSize, c; 17670034214SMatthew G. Knepley PetscErrorCode ierr; 17770034214SMatthew G. Knepley 17870034214SMatthew G. Knepley PetscFunctionBeginHot; 17970034214SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 18070034214SMatthew G. Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 1814b6b44bdSMatthew G. Knepley for (c = 0; c <= coneSize; ++c) { 1824b6b44bdSMatthew G. Knepley const PetscInt point = !c ? p : cone[c-1]; 18370034214SMatthew G. Knepley const PetscInt *support = NULL; 18470034214SMatthew G. Knepley PetscInt supportSize, s, q; 18570034214SMatthew G. Knepley 1864b6b44bdSMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 1874b6b44bdSMatthew G. Knepley ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 18870034214SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 18970034214SMatthew G. Knepley for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) { 19070034214SMatthew G. Knepley if (support[s] == adj[q]) break; 19170034214SMatthew G. Knepley } 19270034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 19370034214SMatthew G. Knepley } 19470034214SMatthew G. Knepley } 19570034214SMatthew G. Knepley *adjSize = numAdj; 19670034214SMatthew G. Knepley PetscFunctionReturn(0); 19770034214SMatthew G. Knepley } 19870034214SMatthew G. Knepley 19970034214SMatthew G. Knepley #undef __FUNCT__ 20070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal" 20170034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 20270034214SMatthew G. Knepley { 20370034214SMatthew G. Knepley const PetscInt *support = NULL; 20470034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, supportSize, s; 20570034214SMatthew G. Knepley PetscErrorCode ierr; 20670034214SMatthew G. Knepley 20770034214SMatthew G. Knepley PetscFunctionBeginHot; 20870034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 20970034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 2104b6b44bdSMatthew G. Knepley for (s = 0; s <= supportSize; ++s) { 2114b6b44bdSMatthew G. Knepley const PetscInt point = !s ? p : support[s-1]; 21270034214SMatthew G. Knepley const PetscInt *cone = NULL; 21370034214SMatthew G. Knepley PetscInt coneSize, c, q; 21470034214SMatthew G. Knepley 2154b6b44bdSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 2164b6b44bdSMatthew G. Knepley ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 21770034214SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 21870034214SMatthew G. Knepley for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) { 21970034214SMatthew G. Knepley if (cone[c] == adj[q]) break; 22070034214SMatthew G. Knepley } 22170034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 22270034214SMatthew G. Knepley } 22370034214SMatthew G. Knepley } 22470034214SMatthew G. Knepley *adjSize = numAdj; 22570034214SMatthew G. Knepley PetscFunctionReturn(0); 22670034214SMatthew G. Knepley } 22770034214SMatthew G. Knepley 22870034214SMatthew G. Knepley #undef __FUNCT__ 22970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal" 23070034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[]) 23170034214SMatthew G. Knepley { 23270034214SMatthew G. Knepley PetscInt *star = NULL; 23370034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, starSize, s; 23470034214SMatthew G. Knepley PetscErrorCode ierr; 23570034214SMatthew G. Knepley 23670034214SMatthew G. Knepley PetscFunctionBeginHot; 23770034214SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 23870034214SMatthew G. Knepley for (s = 0; s < starSize*2; s += 2) { 23970034214SMatthew G. Knepley const PetscInt *closure = NULL; 24070034214SMatthew G. Knepley PetscInt closureSize, c, q; 24170034214SMatthew G. Knepley 24270034214SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 24370034214SMatthew G. Knepley for (c = 0; c < closureSize*2; c += 2) { 24470034214SMatthew G. Knepley for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) { 24570034214SMatthew G. Knepley if (closure[c] == adj[q]) break; 24670034214SMatthew G. Knepley } 24770034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 24870034214SMatthew G. Knepley } 24970034214SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 25070034214SMatthew G. Knepley } 25170034214SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 25270034214SMatthew G. Knepley *adjSize = numAdj; 25370034214SMatthew G. Knepley PetscFunctionReturn(0); 25470034214SMatthew G. Knepley } 25570034214SMatthew G. Knepley 25670034214SMatthew G. Knepley #undef __FUNCT__ 25770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal" 2585b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[]) 25970034214SMatthew G. Knepley { 26079bad331SMatthew G. Knepley static PetscInt asiz = 0; 2618b0b4c70SToby Isaac PetscInt maxAnchors = 1; 2628b0b4c70SToby Isaac PetscInt aStart = -1, aEnd = -1; 2638b0b4c70SToby Isaac PetscInt maxAdjSize; 2648b0b4c70SToby Isaac PetscSection aSec = NULL; 2658b0b4c70SToby Isaac IS aIS = NULL; 2668b0b4c70SToby Isaac const PetscInt *anchors; 26770034214SMatthew G. Knepley PetscErrorCode ierr; 26870034214SMatthew G. Knepley 26970034214SMatthew G. Knepley PetscFunctionBeginHot; 2705b317d89SToby Isaac if (useAnchors) { 271a17985deSToby Isaac ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr); 2728b0b4c70SToby Isaac if (aSec) { 2738b0b4c70SToby Isaac ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr); 27424c766afSToby Isaac maxAnchors = PetscMax(1,maxAnchors); 2758b0b4c70SToby Isaac ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr); 2768b0b4c70SToby Isaac ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 2778b0b4c70SToby Isaac } 2788b0b4c70SToby Isaac } 27979bad331SMatthew G. Knepley if (!*adj) { 28024c766afSToby Isaac PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd; 28179bad331SMatthew G. Knepley 28224c766afSToby Isaac ierr = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr); 28379bad331SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 28424c766afSToby Isaac ierr = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr); 28524c766afSToby Isaac coneSeries = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1; 28624c766afSToby Isaac supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1; 28724c766afSToby Isaac asiz = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries); 2888b0b4c70SToby Isaac asiz *= maxAnchors; 28924c766afSToby Isaac asiz = PetscMin(asiz,pEnd-pStart); 29079bad331SMatthew G. Knepley ierr = PetscMalloc1(asiz,adj);CHKERRQ(ierr); 29179bad331SMatthew G. Knepley } 29279bad331SMatthew G. Knepley if (*adjSize < 0) *adjSize = asiz; 2938b0b4c70SToby Isaac maxAdjSize = *adjSize; 29470034214SMatthew G. Knepley if (useTransitiveClosure) { 29579bad331SMatthew G. Knepley ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr); 29670034214SMatthew G. Knepley } else if (useCone) { 29779bad331SMatthew G. Knepley ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr); 29870034214SMatthew G. Knepley } else { 29979bad331SMatthew G. Knepley ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr); 30070034214SMatthew G. Knepley } 3015b317d89SToby Isaac if (useAnchors && aSec) { 3028b0b4c70SToby Isaac PetscInt origSize = *adjSize; 3038b0b4c70SToby Isaac PetscInt numAdj = origSize; 3048b0b4c70SToby Isaac PetscInt i = 0, j; 3058b0b4c70SToby Isaac PetscInt *orig = *adj; 3068b0b4c70SToby Isaac 3078b0b4c70SToby Isaac while (i < origSize) { 3088b0b4c70SToby Isaac PetscInt p = orig[i]; 3098b0b4c70SToby Isaac PetscInt aDof = 0; 3108b0b4c70SToby Isaac 3118b0b4c70SToby Isaac if (p >= aStart && p < aEnd) { 3128b0b4c70SToby Isaac ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 3138b0b4c70SToby Isaac } 3148b0b4c70SToby Isaac if (aDof) { 3158b0b4c70SToby Isaac PetscInt aOff; 3168b0b4c70SToby Isaac PetscInt s, q; 3178b0b4c70SToby Isaac 3188b0b4c70SToby Isaac for (j = i + 1; j < numAdj; j++) { 3198b0b4c70SToby Isaac orig[j - 1] = orig[j]; 3208b0b4c70SToby Isaac } 3218b0b4c70SToby Isaac origSize--; 3228b0b4c70SToby Isaac numAdj--; 3238b0b4c70SToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 3248b0b4c70SToby Isaac for (s = 0; s < aDof; ++s) { 3258b0b4c70SToby Isaac for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) { 3268b0b4c70SToby Isaac if (anchors[aOff+s] == orig[q]) break; 3278b0b4c70SToby Isaac } 3288b0b4c70SToby Isaac if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 3298b0b4c70SToby Isaac } 3308b0b4c70SToby Isaac } 3318b0b4c70SToby Isaac else { 3328b0b4c70SToby Isaac i++; 3338b0b4c70SToby Isaac } 3348b0b4c70SToby Isaac } 3358b0b4c70SToby Isaac *adjSize = numAdj; 3368b0b4c70SToby Isaac ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 3378b0b4c70SToby Isaac } 33870034214SMatthew G. Knepley PetscFunctionReturn(0); 33970034214SMatthew G. Knepley } 34070034214SMatthew G. Knepley 34170034214SMatthew G. Knepley #undef __FUNCT__ 34270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency" 34370034214SMatthew G. Knepley /*@ 34470034214SMatthew G. Knepley DMPlexGetAdjacency - Return all points adjacent to the given point 34570034214SMatthew G. Knepley 34670034214SMatthew G. Knepley Input Parameters: 34770034214SMatthew G. Knepley + dm - The DM object 34870034214SMatthew G. Knepley . p - The point 34970034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE 35070034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize 35170034214SMatthew G. Knepley 35270034214SMatthew G. Knepley Output Parameters: 35370034214SMatthew G. Knepley + adjSize - The number of adjacent points 35470034214SMatthew G. Knepley - adj - The adjacent points 35570034214SMatthew G. Knepley 35670034214SMatthew G. Knepley Level: advanced 35770034214SMatthew G. Knepley 35870034214SMatthew G. Knepley Notes: The user must PetscFree the adj array if it was not passed in. 35970034214SMatthew G. Knepley 36070034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator() 36170034214SMatthew G. Knepley @*/ 36270034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[]) 36370034214SMatthew G. Knepley { 36470034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 36570034214SMatthew G. Knepley PetscErrorCode ierr; 36670034214SMatthew G. Knepley 36770034214SMatthew G. Knepley PetscFunctionBeginHot; 36870034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 36970034214SMatthew G. Knepley PetscValidPointer(adjSize,3); 37070034214SMatthew G. Knepley PetscValidPointer(adj,4); 3715b317d89SToby Isaac ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr); 37270034214SMatthew G. Knepley PetscFunctionReturn(0); 37370034214SMatthew G. Knepley } 374b0a623aaSMatthew G. Knepley #undef __FUNCT__ 375b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF" 376b0a623aaSMatthew G. Knepley /*@ 377b0a623aaSMatthew G. Knepley DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity 378b0a623aaSMatthew G. Knepley 379b0a623aaSMatthew G. Knepley Collective on DM 380b0a623aaSMatthew G. Knepley 381b0a623aaSMatthew G. Knepley Input Parameters: 382b0a623aaSMatthew G. Knepley + dm - The DM 383b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity 384b0a623aaSMatthew G. Knepley 385b0a623aaSMatthew G. Knepley Output Parameters: 386b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL 387b0a623aaSMatthew G. Knepley - sfProcess - An SF encoding the two-sided process connectivity, or NULL 388b0a623aaSMatthew G. Knepley 389b0a623aaSMatthew G. Knepley Level: developer 390b0a623aaSMatthew G. Knepley 391b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF() 392b0a623aaSMatthew G. Knepley @*/ 393b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess) 394b0a623aaSMatthew G. Knepley { 395b0a623aaSMatthew G. Knepley const PetscSFNode *remotePoints; 396b0a623aaSMatthew G. Knepley PetscInt *localPointsNew; 397b0a623aaSMatthew G. Knepley PetscSFNode *remotePointsNew; 398b0a623aaSMatthew G. Knepley const PetscInt *nranks; 399b0a623aaSMatthew G. Knepley PetscInt *ranksNew; 400b0a623aaSMatthew G. Knepley PetscBT neighbors; 401b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, numLeaves, l, numNeighbors, n; 402b0a623aaSMatthew G. Knepley PetscMPIInt numProcs, proc, rank; 403b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 404b0a623aaSMatthew G. Knepley 405b0a623aaSMatthew G. Knepley PetscFunctionBegin; 406b0a623aaSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 407b0a623aaSMatthew G. Knepley PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 408b0a623aaSMatthew G. Knepley if (processRanks) {PetscValidPointer(processRanks, 3);} 409b0a623aaSMatthew G. Knepley if (sfProcess) {PetscValidPointer(sfProcess, 4);} 410b0a623aaSMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 411b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 412b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr); 413b0a623aaSMatthew G. Knepley ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr); 414b0a623aaSMatthew G. Knepley ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr); 415b0a623aaSMatthew G. Knepley /* Compute root-to-leaf process connectivity */ 416b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr); 417b0a623aaSMatthew G. Knepley ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr); 418b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 419b0a623aaSMatthew G. Knepley PetscInt ndof, noff, n; 420b0a623aaSMatthew G. Knepley 421b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr); 422b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr); 423b0a623aaSMatthew G. Knepley for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);} 424b0a623aaSMatthew G. Knepley } 425b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr); 426b0a623aaSMatthew G. Knepley /* Compute leaf-to-neighbor process connectivity */ 427b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr); 428b0a623aaSMatthew G. Knepley ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr); 429b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 430b0a623aaSMatthew G. Knepley PetscInt ndof, noff, n; 431b0a623aaSMatthew G. Knepley 432b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr); 433b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr); 434b0a623aaSMatthew G. Knepley for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);} 435b0a623aaSMatthew G. Knepley } 436b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr); 437b0a623aaSMatthew G. Knepley /* Compute leaf-to-root process connectivity */ 438b0a623aaSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);} 439b0a623aaSMatthew G. Knepley /* Calculate edges */ 440b0a623aaSMatthew G. Knepley PetscBTClear(neighbors, rank); 441b0a623aaSMatthew G. Knepley for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;} 442b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr); 443b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr); 444b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr); 445b0a623aaSMatthew G. Knepley for(proc = 0, n = 0; proc < numProcs; ++proc) { 446b0a623aaSMatthew G. Knepley if (PetscBTLookup(neighbors, proc)) { 447b0a623aaSMatthew G. Knepley ranksNew[n] = proc; 448b0a623aaSMatthew G. Knepley localPointsNew[n] = proc; 449b0a623aaSMatthew G. Knepley remotePointsNew[n].index = rank; 450b0a623aaSMatthew G. Knepley remotePointsNew[n].rank = proc; 451b0a623aaSMatthew G. Knepley ++n; 452b0a623aaSMatthew G. Knepley } 453b0a623aaSMatthew G. Knepley } 454b0a623aaSMatthew G. Knepley ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr); 455b0a623aaSMatthew G. Knepley if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);} 456b0a623aaSMatthew G. Knepley else {ierr = PetscFree(ranksNew);CHKERRQ(ierr);} 457b0a623aaSMatthew G. Knepley if (sfProcess) { 458b0a623aaSMatthew G. Knepley ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr); 459b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr); 460b0a623aaSMatthew G. Knepley ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr); 461b0a623aaSMatthew G. Knepley ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 462b0a623aaSMatthew G. Knepley } 463b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 464b0a623aaSMatthew G. Knepley } 465b0a623aaSMatthew G. Knepley 466b0a623aaSMatthew G. Knepley #undef __FUNCT__ 467b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership" 468b0a623aaSMatthew G. Knepley /*@ 469b0a623aaSMatthew G. Knepley DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF. 470b0a623aaSMatthew G. Knepley 471b0a623aaSMatthew G. Knepley Collective on DM 472b0a623aaSMatthew G. Knepley 473b0a623aaSMatthew G. Knepley Input Parameter: 474b0a623aaSMatthew G. Knepley . dm - The DM 475b0a623aaSMatthew G. Knepley 476b0a623aaSMatthew G. Knepley Output Parameters: 477b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point 478b0a623aaSMatthew G. Knepley . rootrank - The rank of each edge into the root point 479b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point 480b0a623aaSMatthew G. Knepley - leafrank - The rank of each process sharing a leaf point 481b0a623aaSMatthew G. Knepley 482b0a623aaSMatthew G. Knepley Level: developer 483b0a623aaSMatthew G. Knepley 484b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap() 485b0a623aaSMatthew G. Knepley @*/ 486b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank) 487b0a623aaSMatthew G. Knepley { 488b0a623aaSMatthew G. Knepley MPI_Comm comm; 489b0a623aaSMatthew G. Knepley PetscSF sfPoint; 490b0a623aaSMatthew G. Knepley const PetscInt *rootdegree; 491b0a623aaSMatthew G. Knepley PetscInt *myrank, *remoterank; 492b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, nedges; 493b0a623aaSMatthew G. Knepley PetscMPIInt rank; 494b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 495b0a623aaSMatthew G. Knepley 496b0a623aaSMatthew G. Knepley PetscFunctionBegin; 497b0a623aaSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 498b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 499b0a623aaSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 500b0a623aaSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 501b0a623aaSMatthew G. Knepley /* Compute number of leaves for each root */ 502b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr); 503b0a623aaSMatthew G. Knepley ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr); 504b0a623aaSMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr); 505b0a623aaSMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr); 506b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);} 507b0a623aaSMatthew G. Knepley ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 508b0a623aaSMatthew G. Knepley /* Gather rank of each leaf to root */ 509b0a623aaSMatthew G. Knepley ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr); 510b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr); 511b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(nedges, &remoterank);CHKERRQ(ierr); 512b0a623aaSMatthew G. Knepley for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank; 513b0a623aaSMatthew G. Knepley ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 514b0a623aaSMatthew G. Knepley ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 515b0a623aaSMatthew G. Knepley ierr = PetscFree(myrank);CHKERRQ(ierr); 516b0a623aaSMatthew G. Knepley ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr); 517b0a623aaSMatthew G. Knepley /* Distribute remote ranks to leaves */ 518b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr); 519b0a623aaSMatthew G. Knepley ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr); 520b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 521b0a623aaSMatthew G. Knepley } 522b0a623aaSMatthew G. Knepley 523b0a623aaSMatthew G. Knepley #undef __FUNCT__ 524b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap" 525278397a0SMatthew G. Knepley /*@C 526b0a623aaSMatthew G. Knepley DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF. 527b0a623aaSMatthew G. Knepley 528b0a623aaSMatthew G. Knepley Collective on DM 529b0a623aaSMatthew G. Knepley 530b0a623aaSMatthew G. Knepley Input Parameters: 531b0a623aaSMatthew G. Knepley + dm - The DM 532b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point 533b0a623aaSMatthew G. Knepley . rootrank - The rank of each edge into the root point 534b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point 535b0a623aaSMatthew G. Knepley - leafrank - The rank of each process sharing a leaf point 536b0a623aaSMatthew G. Knepley 537b0a623aaSMatthew G. Knepley Output Parameters: 538*1fd9873aSMichael Lange + overlapSF - The star forest mapping remote overlap points into a contiguous leaf space 539b0a623aaSMatthew G. Knepley 540b0a623aaSMatthew G. Knepley Level: developer 541b0a623aaSMatthew G. Knepley 542*1fd9873aSMichael Lange .seealso: DMPlexDistributeOwnership(), DMPlexDistribute() 543b0a623aaSMatthew G. Knepley @*/ 544e540f424SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, PetscSF *overlapSF) 545b0a623aaSMatthew G. Knepley { 546e540f424SMichael Lange MPI_Comm comm; 547b0a623aaSMatthew G. Knepley DMLabel ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */ 548b0a623aaSMatthew G. Knepley PetscSF sfPoint, sfProc; 549e540f424SMichael Lange DMLabel ovLeafLabel; 550b0a623aaSMatthew G. Knepley const PetscSFNode *remote; 551b0a623aaSMatthew G. Knepley const PetscInt *local; 552*1fd9873aSMichael Lange const PetscInt *nrank, *rrank; 553b0a623aaSMatthew G. Knepley PetscInt *adj = NULL; 554*1fd9873aSMichael Lange PetscInt pStart, pEnd, p, sStart, sEnd, nleaves, l; 555b0a623aaSMatthew G. Knepley PetscMPIInt rank, numProcs; 55626a7d390SMatthew G. Knepley PetscBool useCone, useClosure, flg; 557b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 558b0a623aaSMatthew G. Knepley 559b0a623aaSMatthew G. Knepley PetscFunctionBegin; 560e540f424SMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 561e540f424SMichael Lange ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 562e540f424SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 563b0a623aaSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 564b0a623aaSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 565b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr); 566b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 567b0a623aaSMatthew G. Knepley ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr); 568b0a623aaSMatthew G. Knepley /* Handle leaves: shared with the root point */ 569b0a623aaSMatthew G. Knepley for (l = 0; l < nleaves; ++l) { 570b0a623aaSMatthew G. Knepley PetscInt adjSize = PETSC_DETERMINE, a; 571b0a623aaSMatthew G. Knepley 572b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr); 573b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);} 574b0a623aaSMatthew G. Knepley } 575b0a623aaSMatthew G. Knepley ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr); 576b0a623aaSMatthew G. Knepley ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr); 577b0a623aaSMatthew G. Knepley /* Handle roots */ 578b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 579b0a623aaSMatthew G. Knepley PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a; 580b0a623aaSMatthew G. Knepley 581b0a623aaSMatthew G. Knepley if ((p >= sStart) && (p < sEnd)) { 582b0a623aaSMatthew G. Knepley /* Some leaves share a root with other leaves on different processes */ 583b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr); 584b0a623aaSMatthew G. Knepley if (neighbors) { 585b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr); 586b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 587b0a623aaSMatthew G. Knepley for (n = 0; n < neighbors; ++n) { 588b0a623aaSMatthew G. Knepley const PetscInt remoteRank = nrank[noff+n]; 589b0a623aaSMatthew G. Knepley 590b0a623aaSMatthew G. Knepley if (remoteRank == rank) continue; 591b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 592b0a623aaSMatthew G. Knepley } 593b0a623aaSMatthew G. Knepley } 594b0a623aaSMatthew G. Knepley } 595b0a623aaSMatthew G. Knepley /* Roots are shared with leaves */ 596b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr); 597b0a623aaSMatthew G. Knepley if (!neighbors) continue; 598b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr); 599b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 600b0a623aaSMatthew G. Knepley for (n = 0; n < neighbors; ++n) { 601b0a623aaSMatthew G. Knepley const PetscInt remoteRank = rrank[noff+n]; 602b0a623aaSMatthew G. Knepley 603b0a623aaSMatthew G. Knepley if (remoteRank == rank) continue; 604b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 605b0a623aaSMatthew G. Knepley } 606b0a623aaSMatthew G. Knepley } 607b0a623aaSMatthew G. Knepley ierr = PetscFree(adj);CHKERRQ(ierr); 608b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr); 609b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr); 61026a7d390SMatthew G. Knepley /* We require the closure in the overlap */ 61126a7d390SMatthew G. Knepley ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr); 61226a7d390SMatthew G. Knepley ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr); 61326a7d390SMatthew G. Knepley if (useCone || !useClosure) { 61426a7d390SMatthew G. Knepley IS rankIS, pointIS; 61526a7d390SMatthew G. Knepley const PetscInt *ranks, *points; 61626a7d390SMatthew G. Knepley PetscInt numRanks, numPoints, r, p; 61726a7d390SMatthew G. Knepley 61826a7d390SMatthew G. Knepley ierr = DMLabelGetValueIS(ovAdjByRank, &rankIS);CHKERRQ(ierr); 61926a7d390SMatthew G. Knepley ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 62026a7d390SMatthew G. Knepley ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 62126a7d390SMatthew G. Knepley for (r = 0; r < numRanks; ++r) { 62226a7d390SMatthew G. Knepley const PetscInt rank = ranks[r]; 62326a7d390SMatthew G. Knepley 62426a7d390SMatthew G. Knepley ierr = DMLabelGetStratumIS(ovAdjByRank, rank, &pointIS);CHKERRQ(ierr); 62526a7d390SMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 62626a7d390SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 62726a7d390SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 62826a7d390SMatthew G. Knepley PetscInt *closure = NULL, closureSize, c; 62926a7d390SMatthew G. Knepley 63026a7d390SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 63126a7d390SMatthew G. Knepley for (c = 0; c < closureSize*2; c += 2) {ierr = DMLabelSetValue(ovAdjByRank, closure[c], rank);CHKERRQ(ierr);} 63226a7d390SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 63326a7d390SMatthew G. Knepley } 63426a7d390SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 63526a7d390SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 63626a7d390SMatthew G. Knepley } 63726a7d390SMatthew G. Knepley ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 63826a7d390SMatthew G. Knepley ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 63926a7d390SMatthew G. Knepley } 640e540f424SMichael Lange ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr); 641e540f424SMichael Lange if (flg) { 642b0a623aaSMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 643b0a623aaSMatthew G. Knepley ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 644b0a623aaSMatthew G. Knepley } 645*1fd9873aSMichael Lange /* Make process SF and invert sender to receiver label */ 646b0a623aaSMatthew G. Knepley ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr); 647e540f424SMichael Lange ierr = DMLabelCreate("ovLeafs", &ovLeafLabel);CHKERRQ(ierr); 648*1fd9873aSMichael Lange ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, ovLeafLabel);CHKERRQ(ierr); 649e540f424SMichael Lange /* Don't import points from yourself */ 650*1fd9873aSMichael Lange ierr = DMLabelClearStratum(ovLeafLabel, rank);CHKERRQ(ierr); 651e540f424SMichael Lange /* Don't import points already in the pointSF */ 652e540f424SMichael Lange for (l = 0; l < nleaves; ++l) { 653e540f424SMichael Lange ierr = DMLabelClearValue(ovLeafLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr); 654e540f424SMichael Lange } 655aa3148a8SMichael Lange /* Convert receiver label to SF */ 656aa3148a8SMichael Lange ierr = DMPlexPartitionLabelCreateSF(dm, ovLeafLabel, overlapSF);CHKERRQ(ierr); 657e540f424SMichael Lange ierr = PetscObjectSetName((PetscObject) *overlapSF, "Overlap SF");CHKERRQ(ierr); 658e540f424SMichael Lange ierr = PetscSFSetFromOptions(*overlapSF);CHKERRQ(ierr); 659e540f424SMichael Lange if (flg) { 660e540f424SMichael Lange ierr = PetscPrintf(comm, "Overlap SF\n");CHKERRQ(ierr); 661e540f424SMichael Lange ierr = PetscSFView(*overlapSF, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 662e540f424SMichael Lange } 663e540f424SMichael Lange /* Clean up */ 664*1fd9873aSMichael Lange ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr); 665aa3148a8SMichael Lange ierr = DMLabelDestroy(&ovLeafLabel);CHKERRQ(ierr); 666b0a623aaSMatthew G. Knepley ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr); 667b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 668b0a623aaSMatthew G. Knepley } 66970034214SMatthew G. Knepley 67070034214SMatthew G. Knepley #undef __FUNCT__ 67146f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF" 67246f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF) 67346f9b1c3SMichael Lange { 67446f9b1c3SMichael Lange MPI_Comm comm; 67546f9b1c3SMichael Lange PetscMPIInt rank, numProcs; 67646f9b1c3SMichael Lange PetscInt d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints; 67746f9b1c3SMichael Lange PetscInt *pointDepths, *remoteDepths, *ilocal; 67846f9b1c3SMichael Lange PetscInt *depthRecv, *depthShift, *depthIdx; 67946f9b1c3SMichael Lange PetscSFNode *iremote; 68046f9b1c3SMichael Lange PetscSF pointSF; 68146f9b1c3SMichael Lange const PetscInt *sharedLocal; 68246f9b1c3SMichael Lange const PetscSFNode *overlapRemote, *sharedRemote; 68346f9b1c3SMichael Lange PetscErrorCode ierr; 68446f9b1c3SMichael Lange 68546f9b1c3SMichael Lange PetscFunctionBegin; 68646f9b1c3SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 68746f9b1c3SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 68846f9b1c3SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 68946f9b1c3SMichael Lange ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 69046f9b1c3SMichael Lange ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 69146f9b1c3SMichael Lange 69246f9b1c3SMichael Lange /* Before building the migration SF we need to know the new stratum offsets */ 69346f9b1c3SMichael Lange ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr); 69446f9b1c3SMichael Lange ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr); 69546f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 69646f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 69746f9b1c3SMichael Lange for (p=pStart; p<pEnd; p++) pointDepths[p] = d; 69846f9b1c3SMichael Lange } 69946f9b1c3SMichael Lange for (p=0; p<nleaves; p++) remoteDepths[p] = -1; 70046f9b1c3SMichael Lange ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 70146f9b1c3SMichael Lange ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 70246f9b1c3SMichael Lange 70346f9b1c3SMichael Lange /* Count recevied points in each stratum and compute the internal strata shift */ 70446f9b1c3SMichael Lange ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr); 70546f9b1c3SMichael Lange for (d=0; d<dim+1; d++) depthRecv[d]=0; 70646f9b1c3SMichael Lange for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++; 70746f9b1c3SMichael Lange depthShift[dim] = 0; 70846f9b1c3SMichael Lange for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim]; 70946f9b1c3SMichael Lange for (d=1; d<dim; d++) depthShift[d] += depthRecv[0]; 71046f9b1c3SMichael Lange for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1]; 71146f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 71246f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 71346f9b1c3SMichael Lange depthIdx[d] = pStart + depthShift[d]; 71446f9b1c3SMichael Lange } 71546f9b1c3SMichael Lange 71646f9b1c3SMichael Lange /* Form the overlap SF build an SF that describes the full overlap migration SF */ 71746f9b1c3SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 71846f9b1c3SMichael Lange newLeaves = pEnd - pStart + nleaves; 71909b7985cSMatthew G. Knepley ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr); 72009b7985cSMatthew G. Knepley ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr); 72146f9b1c3SMichael Lange /* First map local points to themselves */ 72246f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 72346f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 72446f9b1c3SMichael Lange for (p=pStart; p<pEnd; p++) { 72546f9b1c3SMichael Lange point = p + depthShift[d]; 72646f9b1c3SMichael Lange ilocal[point] = point; 72746f9b1c3SMichael Lange iremote[point].index = p; 72846f9b1c3SMichael Lange iremote[point].rank = rank; 72946f9b1c3SMichael Lange depthIdx[d]++; 73046f9b1c3SMichael Lange } 73146f9b1c3SMichael Lange } 73246f9b1c3SMichael Lange 73346f9b1c3SMichael Lange /* Add in the remote roots for currently shared points */ 73446f9b1c3SMichael Lange ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr); 73546f9b1c3SMichael Lange ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr); 73646f9b1c3SMichael Lange for (d=0; d<dim+1; d++) { 73746f9b1c3SMichael Lange ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 73846f9b1c3SMichael Lange for (p=0; p<numSharedPoints; p++) { 73946f9b1c3SMichael Lange if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) { 74046f9b1c3SMichael Lange point = sharedLocal[p] + depthShift[d]; 74146f9b1c3SMichael Lange iremote[point].index = sharedRemote[p].index; 74246f9b1c3SMichael Lange iremote[point].rank = sharedRemote[p].rank; 74346f9b1c3SMichael Lange } 74446f9b1c3SMichael Lange } 74546f9b1c3SMichael Lange } 74646f9b1c3SMichael Lange 74746f9b1c3SMichael Lange /* Now add the incoming overlap points */ 74846f9b1c3SMichael Lange for (p=0; p<nleaves; p++) { 74946f9b1c3SMichael Lange point = depthIdx[remoteDepths[p]]; 75046f9b1c3SMichael Lange ilocal[point] = point; 75146f9b1c3SMichael Lange iremote[point].index = overlapRemote[p].index; 75246f9b1c3SMichael Lange iremote[point].rank = overlapRemote[p].rank; 75346f9b1c3SMichael Lange depthIdx[remoteDepths[p]]++; 75446f9b1c3SMichael Lange } 75515fff7beSMatthew G. Knepley ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr); 75646f9b1c3SMichael Lange 75746f9b1c3SMichael Lange ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr); 75846f9b1c3SMichael Lange ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr); 75946f9b1c3SMichael Lange ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr); 76046f9b1c3SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 76146f9b1c3SMichael Lange ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 76246f9b1c3SMichael Lange 76346f9b1c3SMichael Lange ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr); 76470034214SMatthew G. Knepley PetscFunctionReturn(0); 76570034214SMatthew G. Knepley } 76670034214SMatthew G. Knepley 76770034214SMatthew G. Knepley #undef __FUNCT__ 76870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField" 76970034214SMatthew G. Knepley /*@ 77070034214SMatthew G. Knepley DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 77170034214SMatthew G. Knepley 77270034214SMatthew G. Knepley Collective on DM 77370034214SMatthew G. Knepley 77470034214SMatthew G. Knepley Input Parameters: 77570034214SMatthew G. Knepley + dm - The DMPlex object 77670034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 77770034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 77870034214SMatthew G. Knepley - originalVec - The existing data 77970034214SMatthew G. Knepley 78070034214SMatthew G. Knepley Output Parameters: 78170034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 78270034214SMatthew G. Knepley - newVec - The new data 78370034214SMatthew G. Knepley 78470034214SMatthew G. Knepley Level: developer 78570034214SMatthew G. Knepley 78670034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData() 78770034214SMatthew G. Knepley @*/ 78870034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec) 78970034214SMatthew G. Knepley { 79070034214SMatthew G. Knepley PetscSF fieldSF; 79170034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 79270034214SMatthew G. Knepley PetscScalar *originalValues, *newValues; 79370034214SMatthew G. Knepley PetscErrorCode ierr; 79470034214SMatthew G. Knepley 79570034214SMatthew G. Knepley PetscFunctionBegin; 79670034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 79770034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 79870034214SMatthew G. Knepley 79970034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 80070034214SMatthew G. Knepley ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr); 80170034214SMatthew G. Knepley ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr); 80270034214SMatthew G. Knepley 80370034214SMatthew G. Knepley ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr); 80470034214SMatthew G. Knepley ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr); 80570034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 80670034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 80770034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 80870034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 80970034214SMatthew G. Knepley ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr); 81070034214SMatthew G. Knepley ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr); 81170034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 81270034214SMatthew G. Knepley PetscFunctionReturn(0); 81370034214SMatthew G. Knepley } 81470034214SMatthew G. Knepley 81570034214SMatthew G. Knepley #undef __FUNCT__ 81670034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS" 81770034214SMatthew G. Knepley /*@ 81870034214SMatthew G. Knepley DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 81970034214SMatthew G. Knepley 82070034214SMatthew G. Knepley Collective on DM 82170034214SMatthew G. Knepley 82270034214SMatthew G. Knepley Input Parameters: 82370034214SMatthew G. Knepley + dm - The DMPlex object 82470034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 82570034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 82670034214SMatthew G. Knepley - originalIS - The existing data 82770034214SMatthew G. Knepley 82870034214SMatthew G. Knepley Output Parameters: 82970034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 83070034214SMatthew G. Knepley - newIS - The new data 83170034214SMatthew G. Knepley 83270034214SMatthew G. Knepley Level: developer 83370034214SMatthew G. Knepley 83470034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData() 83570034214SMatthew G. Knepley @*/ 83670034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS) 83770034214SMatthew G. Knepley { 83870034214SMatthew G. Knepley PetscSF fieldSF; 83970034214SMatthew G. Knepley PetscInt *newValues, *remoteOffsets, fieldSize; 84070034214SMatthew G. Knepley const PetscInt *originalValues; 84170034214SMatthew G. Knepley PetscErrorCode ierr; 84270034214SMatthew G. Knepley 84370034214SMatthew G. Knepley PetscFunctionBegin; 84470034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 84570034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 84670034214SMatthew G. Knepley 84770034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 84870034214SMatthew G. Knepley ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr); 84970034214SMatthew G. Knepley 85070034214SMatthew G. Knepley ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr); 85170034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 852aaf8c182SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr); 853aaf8c182SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr); 85470034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 85570034214SMatthew G. Knepley ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr); 85670034214SMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr); 85770034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 85870034214SMatthew G. Knepley PetscFunctionReturn(0); 85970034214SMatthew G. Knepley } 86070034214SMatthew G. Knepley 86170034214SMatthew G. Knepley #undef __FUNCT__ 86270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData" 86370034214SMatthew G. Knepley /*@ 86470034214SMatthew G. Knepley DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 86570034214SMatthew G. Knepley 86670034214SMatthew G. Knepley Collective on DM 86770034214SMatthew G. Knepley 86870034214SMatthew G. Knepley Input Parameters: 86970034214SMatthew G. Knepley + dm - The DMPlex object 87070034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 87170034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 87270034214SMatthew G. Knepley . datatype - The type of data 87370034214SMatthew G. Knepley - originalData - The existing data 87470034214SMatthew G. Knepley 87570034214SMatthew G. Knepley Output Parameters: 87670034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout 87770034214SMatthew G. Knepley - newData - The new data 87870034214SMatthew G. Knepley 87970034214SMatthew G. Knepley Level: developer 88070034214SMatthew G. Knepley 88170034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField() 88270034214SMatthew G. Knepley @*/ 88370034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData) 88470034214SMatthew G. Knepley { 88570034214SMatthew G. Knepley PetscSF fieldSF; 88670034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 88770034214SMatthew G. Knepley PetscMPIInt dataSize; 88870034214SMatthew G. Knepley PetscErrorCode ierr; 88970034214SMatthew G. Knepley 89070034214SMatthew G. Knepley PetscFunctionBegin; 89170034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 89270034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 89370034214SMatthew G. Knepley 89470034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 89570034214SMatthew G. Knepley ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr); 89670034214SMatthew G. Knepley ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr); 89770034214SMatthew G. Knepley 89870034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 89970034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 90070034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 90170034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 90270034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 90370034214SMatthew G. Knepley PetscFunctionReturn(0); 90470034214SMatthew G. Knepley } 90570034214SMatthew G. Knepley 90670034214SMatthew G. Knepley #undef __FUNCT__ 907cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones" 908cc71bff1SMichael Lange PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 909cc71bff1SMichael Lange { 910cc71bff1SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 911cc71bff1SMichael Lange MPI_Comm comm; 912cc71bff1SMichael Lange PetscSF coneSF; 913cc71bff1SMichael Lange PetscSection originalConeSection, newConeSection; 914cc71bff1SMichael Lange PetscInt *remoteOffsets, *cones, *newCones, newConesSize; 915cc71bff1SMichael Lange PetscBool flg; 916cc71bff1SMichael Lange PetscErrorCode ierr; 917cc71bff1SMichael Lange 918cc71bff1SMichael Lange PetscFunctionBegin; 919cc71bff1SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 920cc71bff1SMichael Lange PetscValidPointer(dmParallel,4); 921cc71bff1SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 922cc71bff1SMichael Lange 923cc71bff1SMichael Lange /* Distribute cone section */ 924cc71bff1SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 925cc71bff1SMichael Lange ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); 926cc71bff1SMichael Lange ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr); 927cc71bff1SMichael Lange ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); 928cc71bff1SMichael Lange ierr = DMSetUp(dmParallel);CHKERRQ(ierr); 929cc71bff1SMichael Lange { 930cc71bff1SMichael Lange PetscInt pStart, pEnd, p; 931cc71bff1SMichael Lange 932cc71bff1SMichael Lange ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr); 933cc71bff1SMichael Lange for (p = pStart; p < pEnd; ++p) { 934cc71bff1SMichael Lange PetscInt coneSize; 935cc71bff1SMichael Lange ierr = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr); 936cc71bff1SMichael Lange pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize); 937cc71bff1SMichael Lange } 938cc71bff1SMichael Lange } 939cc71bff1SMichael Lange /* Communicate and renumber cones */ 940cc71bff1SMichael Lange ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); 941cc71bff1SMichael Lange ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 942cc71bff1SMichael Lange ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr); 943cc71bff1SMichael Lange ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 944cc71bff1SMichael Lange ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 945cc71bff1SMichael Lange ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); 946cc71bff1SMichael Lange ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr); 9473533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG 9483533c52bSMatthew G. Knepley { 9493533c52bSMatthew G. Knepley PetscInt p; 9503533c52bSMatthew G. Knepley PetscBool valid = PETSC_TRUE; 9513533c52bSMatthew G. Knepley for (p = 0; p < newConesSize; ++p) { 9523533c52bSMatthew G. Knepley if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);} 9533533c52bSMatthew G. Knepley } 9543533c52bSMatthew G. Knepley if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map"); 9553533c52bSMatthew G. Knepley } 9563533c52bSMatthew G. Knepley #endif 957cc71bff1SMichael Lange ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr); 958cc71bff1SMichael Lange if (flg) { 959cc71bff1SMichael Lange ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); 960cc71bff1SMichael Lange ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 961cc71bff1SMichael Lange ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); 962cc71bff1SMichael Lange ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 963cc71bff1SMichael Lange ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr); 964cc71bff1SMichael Lange } 965cc71bff1SMichael Lange ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr); 966cc71bff1SMichael Lange ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr); 967cc71bff1SMichael Lange ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 968cc71bff1SMichael Lange ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 969cc71bff1SMichael Lange ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); 970cc71bff1SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 971cc71bff1SMichael Lange /* Create supports and stratify sieve */ 972cc71bff1SMichael Lange { 973cc71bff1SMichael Lange PetscInt pStart, pEnd; 974cc71bff1SMichael Lange 975cc71bff1SMichael Lange ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 976cc71bff1SMichael Lange ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr); 977cc71bff1SMichael Lange } 978cc71bff1SMichael Lange ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr); 979cc71bff1SMichael Lange ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr); 980cc71bff1SMichael Lange PetscFunctionReturn(0); 981cc71bff1SMichael Lange } 982cc71bff1SMichael Lange 983cc71bff1SMichael Lange #undef __FUNCT__ 9840df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates" 9850df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel) 9860df0e737SMichael Lange { 9870df0e737SMichael Lange MPI_Comm comm; 9880df0e737SMichael Lange PetscSection originalCoordSection, newCoordSection; 9890df0e737SMichael Lange Vec originalCoordinates, newCoordinates; 9900df0e737SMichael Lange PetscInt bs; 9910df0e737SMichael Lange const char *name; 9920df0e737SMichael Lange const PetscReal *maxCell, *L; 9930df0e737SMichael Lange PetscErrorCode ierr; 9940df0e737SMichael Lange 9950df0e737SMichael Lange PetscFunctionBegin; 9960df0e737SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9970df0e737SMichael Lange PetscValidPointer(dmParallel, 3); 9980df0e737SMichael Lange 9990df0e737SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 10000df0e737SMichael Lange ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); 10010df0e737SMichael Lange ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr); 10020df0e737SMichael Lange ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr); 10030df0e737SMichael Lange if (originalCoordinates) { 10040df0e737SMichael Lange ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr); 10050df0e737SMichael Lange ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr); 10060df0e737SMichael Lange ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr); 10070df0e737SMichael Lange 10080df0e737SMichael Lange ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr); 10090df0e737SMichael Lange ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr); 10100df0e737SMichael Lange ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr); 10110df0e737SMichael Lange ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr); 10120df0e737SMichael Lange ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 10130df0e737SMichael Lange } 10140df0e737SMichael Lange ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr); 10150df0e737SMichael Lange if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L);CHKERRQ(ierr);} 10160df0e737SMichael Lange PetscFunctionReturn(0); 10170df0e737SMichael Lange } 10180df0e737SMichael Lange 10190df0e737SMichael Lange #undef __FUNCT__ 10200df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels" 1021d995df53SMatthew G. Knepley /* Here we are assuming that process 0 always has everything */ 10222eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel) 10230df0e737SMichael Lange { 10240df0e737SMichael Lange MPI_Comm comm; 10250df0e737SMichael Lange PetscMPIInt rank; 1026d995df53SMatthew G. Knepley PetscInt numLabels, numLocalLabels, l; 1027d995df53SMatthew G. Knepley PetscBool hasLabels = PETSC_FALSE; 10280df0e737SMichael Lange PetscErrorCode ierr; 10290df0e737SMichael Lange 10300df0e737SMichael Lange PetscFunctionBegin; 10310df0e737SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10322eb1fa14SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 10330df0e737SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 10340df0e737SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 10350df0e737SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10360df0e737SMichael Lange 1037d995df53SMatthew G. Knepley /* Everyone must have either the same number of labels, or none */ 1038d995df53SMatthew G. Knepley ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr); 1039d995df53SMatthew G. Knepley numLabels = numLocalLabels; 10400df0e737SMichael Lange ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1041627847f0SMatthew G. Knepley if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE; 1042bdd2d751SMichael Lange for (l = numLabels-1; l >= 0; --l) { 1043bdd2d751SMichael Lange DMLabel label = NULL, labelNew = NULL; 10440df0e737SMichael Lange PetscBool isdepth; 10450df0e737SMichael Lange 1046d995df53SMatthew G. Knepley if (hasLabels) { 1047bdd2d751SMichael Lange ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr); 10480df0e737SMichael Lange /* Skip "depth" because it is recreated */ 1049bdd2d751SMichael Lange ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr); 1050bdd2d751SMichael Lange } 10510df0e737SMichael Lange ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 1052bdd2d751SMichael Lange if (isdepth) continue; 1053bdd2d751SMichael Lange ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr); 1054bdd2d751SMichael Lange ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr); 10550df0e737SMichael Lange } 10560df0e737SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 10570df0e737SMichael Lange PetscFunctionReturn(0); 10580df0e737SMichael Lange } 10590df0e737SMichael Lange 10609734c634SMichael Lange #undef __FUNCT__ 10619734c634SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupHybrid" 10629734c634SMichael Lange PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 10639734c634SMichael Lange { 10649734c634SMichael Lange DM_Plex *mesh = (DM_Plex*) dm->data; 10659734c634SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 10669734c634SMichael Lange MPI_Comm comm; 10679734c634SMichael Lange const PetscInt *gpoints; 10689734c634SMichael Lange PetscInt dim, depth, n, d; 10699734c634SMichael Lange PetscErrorCode ierr; 10709734c634SMichael Lange 10719734c634SMichael Lange PetscFunctionBegin; 10729734c634SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10739734c634SMichael Lange PetscValidPointer(dmParallel, 4); 10749734c634SMichael Lange 10759734c634SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 10769734c634SMichael Lange ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 10779734c634SMichael Lange 10789734c634SMichael Lange /* Setup hybrid structure */ 10799734c634SMichael Lange for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];} 10809734c634SMichael Lange ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr); 10819734c634SMichael Lange ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr); 10829734c634SMichael Lange ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr); 10839734c634SMichael Lange ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 10849734c634SMichael Lange for (d = 0; d <= dim; ++d) { 10859734c634SMichael Lange PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p; 10869734c634SMichael Lange 10879734c634SMichael Lange if (pmax < 0) continue; 10889734c634SMichael Lange ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr); 10899734c634SMichael Lange ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr); 10909734c634SMichael Lange ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr); 10919734c634SMichael Lange for (p = 0; p < n; ++p) { 10929734c634SMichael Lange const PetscInt point = gpoints[p]; 10939734c634SMichael Lange 10949734c634SMichael Lange if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax; 10959734c634SMichael Lange } 10969734c634SMichael Lange if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax; 10979734c634SMichael Lange else pmesh->hybridPointMax[d] = -1; 10989734c634SMichael Lange } 10999734c634SMichael Lange ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr); 11009734c634SMichael Lange PetscFunctionReturn(0); 11019734c634SMichael Lange } 11020df0e737SMichael Lange 1103a6f36705SMichael Lange #undef __FUNCT__ 1104a6f36705SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupTree" 1105a6f36705SMichael Lange PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 1106a6f36705SMichael Lange { 1107a6f36705SMichael Lange MPI_Comm comm; 1108a6f36705SMichael Lange DM refTree; 1109a6f36705SMichael Lange PetscSection origParentSection, newParentSection; 1110a6f36705SMichael Lange PetscInt *origParents, *origChildIDs; 1111a6f36705SMichael Lange PetscBool flg; 1112a6f36705SMichael Lange PetscErrorCode ierr; 1113a6f36705SMichael Lange 1114a6f36705SMichael Lange PetscFunctionBegin; 1115a6f36705SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1116a6f36705SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 4); 1117a6f36705SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1118a6f36705SMichael Lange 1119a6f36705SMichael Lange /* Set up tree */ 1120a6f36705SMichael Lange ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 1121a6f36705SMichael Lange ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr); 1122a6f36705SMichael Lange ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr); 1123a6f36705SMichael Lange if (origParentSection) { 1124a6f36705SMichael Lange PetscInt pStart, pEnd; 1125a6f36705SMichael Lange PetscInt *newParents, *newChildIDs; 1126a6f36705SMichael Lange PetscInt *remoteOffsetsParents, newParentSize; 1127a6f36705SMichael Lange PetscSF parentSF; 1128a6f36705SMichael Lange 1129a6f36705SMichael Lange ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 1130a6f36705SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr); 1131a6f36705SMichael Lange ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr); 1132a6f36705SMichael Lange ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr); 1133a6f36705SMichael Lange ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr); 1134a6f36705SMichael Lange ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr); 1135a6f36705SMichael Lange ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr); 1136a6f36705SMichael Lange ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr); 1137a6f36705SMichael Lange ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr); 1138a6f36705SMichael Lange ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr); 1139a6f36705SMichael Lange ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr); 1140a6f36705SMichael Lange ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr); 1141a6f36705SMichael Lange ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr); 1142a6f36705SMichael Lange if (flg) { 1143a6f36705SMichael Lange ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr); 1144a6f36705SMichael Lange ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1145a6f36705SMichael Lange ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr); 1146a6f36705SMichael Lange ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1147a6f36705SMichael Lange ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr); 1148a6f36705SMichael Lange } 1149a6f36705SMichael Lange ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr); 1150a6f36705SMichael Lange ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr); 1151a6f36705SMichael Lange ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr); 1152a6f36705SMichael Lange ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr); 1153a6f36705SMichael Lange } 1154a6f36705SMichael Lange PetscFunctionReturn(0); 1155a6f36705SMichael Lange } 11560df0e737SMichael Lange 11570df0e737SMichael Lange #undef __FUNCT__ 11584eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF" 11594eca1733SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, PetscSection partSection, IS part, PetscSection origPartSection, IS origPart, DM dmParallel) 11604eca1733SMichael Lange { 11614eca1733SMichael Lange DM_Plex *mesh = (DM_Plex*) dm->data; 11624eca1733SMichael Lange DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 11634eca1733SMichael Lange PetscMPIInt rank, numProcs; 11644eca1733SMichael Lange MPI_Comm comm; 11654eca1733SMichael Lange PetscErrorCode ierr; 11664eca1733SMichael Lange 11674eca1733SMichael Lange PetscFunctionBegin; 11684eca1733SMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11694eca1733SMichael Lange PetscValidPointer(dmParallel,7); 11704eca1733SMichael Lange 11714eca1733SMichael Lange /* Create point SF for parallel mesh */ 11724eca1733SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 11734eca1733SMichael Lange ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 11744eca1733SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 11754eca1733SMichael Lange ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 11764eca1733SMichael Lange { 11774eca1733SMichael Lange const PetscInt *leaves; 11784eca1733SMichael Lange PetscSFNode *remotePoints, *rowners, *lowners; 11794eca1733SMichael Lange PetscInt numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints; 11804eca1733SMichael Lange PetscInt pStart, pEnd; 11814eca1733SMichael Lange 11824eca1733SMichael Lange ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 11834eca1733SMichael Lange ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr); 11844eca1733SMichael Lange ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr); 11854eca1733SMichael Lange for (p=0; p<numRoots; p++) { 11864eca1733SMichael Lange rowners[p].rank = -1; 11874eca1733SMichael Lange rowners[p].index = -1; 11884eca1733SMichael Lange } 11894eca1733SMichael Lange if (origPart) { 11904eca1733SMichael Lange /* Make sure points in the original partition are not assigned to other procs */ 11914eca1733SMichael Lange const PetscInt *origPoints; 11924eca1733SMichael Lange 11934eca1733SMichael Lange ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr); 11944eca1733SMichael Lange for (p = 0; p < numProcs; ++p) { 11954eca1733SMichael Lange PetscInt dof, off, d; 11964eca1733SMichael Lange 11974eca1733SMichael Lange ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr); 11984eca1733SMichael Lange ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr); 11994eca1733SMichael Lange for (d = off; d < off+dof; ++d) { 12004eca1733SMichael Lange rowners[origPoints[d]].rank = p; 12014eca1733SMichael Lange } 12024eca1733SMichael Lange } 12034eca1733SMichael Lange ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr); 12044eca1733SMichael Lange } 12054eca1733SMichael Lange ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 12064eca1733SMichael Lange ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 12074eca1733SMichael Lange for (p = 0; p < numLeaves; ++p) { 12084eca1733SMichael Lange if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */ 12094eca1733SMichael Lange lowners[p].rank = rank; 12104eca1733SMichael Lange lowners[p].index = leaves ? leaves[p] : p; 12114eca1733SMichael Lange } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */ 12124eca1733SMichael Lange lowners[p].rank = -2; 12134eca1733SMichael Lange lowners[p].index = -2; 12144eca1733SMichael Lange } 12154eca1733SMichael Lange } 12164eca1733SMichael Lange for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */ 12174eca1733SMichael Lange rowners[p].rank = -3; 12184eca1733SMichael Lange rowners[p].index = -3; 12194eca1733SMichael Lange } 12204eca1733SMichael Lange ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 12214eca1733SMichael Lange ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 12224eca1733SMichael Lange ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 12234eca1733SMichael Lange ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 12244eca1733SMichael Lange for (p = 0; p < numLeaves; ++p) { 12254eca1733SMichael Lange if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed"); 12264eca1733SMichael Lange if (lowners[p].rank != rank) ++numGhostPoints; 12274eca1733SMichael Lange } 12284eca1733SMichael Lange ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr); 12294eca1733SMichael Lange ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr); 12304eca1733SMichael Lange for (p = 0, gp = 0; p < numLeaves; ++p) { 12314eca1733SMichael Lange if (lowners[p].rank != rank) { 12324eca1733SMichael Lange ghostPoints[gp] = leaves ? leaves[p] : p; 12334eca1733SMichael Lange remotePoints[gp].rank = lowners[p].rank; 12344eca1733SMichael Lange remotePoints[gp].index = lowners[p].index; 12354eca1733SMichael Lange ++gp; 12364eca1733SMichael Lange } 12374eca1733SMichael Lange } 12384eca1733SMichael Lange ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr); 12394eca1733SMichael Lange ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 12404eca1733SMichael Lange ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr); 12414eca1733SMichael Lange } 12424eca1733SMichael Lange pmesh->useCone = mesh->useCone; 12434eca1733SMichael Lange pmesh->useClosure = mesh->useClosure; 12444eca1733SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 12454eca1733SMichael Lange PetscFunctionReturn(0); 12464eca1733SMichael Lange } 12474eca1733SMichael Lange 12484eca1733SMichael Lange #undef __FUNCT__ 124970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute" 125070034214SMatthew G. Knepley /*@C 125170034214SMatthew G. Knepley DMPlexDistribute - Distributes the mesh and any associated sections. 125270034214SMatthew G. Knepley 125370034214SMatthew G. Knepley Not Collective 125470034214SMatthew G. Knepley 125570034214SMatthew G. Knepley Input Parameter: 125670034214SMatthew G. Knepley + dm - The original DMPlex object 125770034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default 125870034214SMatthew G. Knepley 125970034214SMatthew G. Knepley Output Parameter: 126070034214SMatthew G. Knepley + sf - The PetscSF used for point distribution 126170034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL 126270034214SMatthew G. Knepley 126370034214SMatthew G. Knepley Note: If the mesh was not distributed, the return value is NULL. 126470034214SMatthew G. Knepley 126570034214SMatthew G. Knepley The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 126670034214SMatthew G. Knepley DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 126770034214SMatthew G. Knepley representation on the mesh. 126870034214SMatthew G. Knepley 126970034214SMatthew G. Knepley Level: intermediate 127070034214SMatthew G. Knepley 127170034214SMatthew G. Knepley .keywords: mesh, elements 127270034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 127370034214SMatthew G. Knepley @*/ 127480cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel) 127570034214SMatthew G. Knepley { 127670034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data, *pmesh; 127770034214SMatthew G. Knepley MPI_Comm comm; 127843331d4aSMichael Lange PetscInt dim, numRemoteRanks, nroots, nleaves; 1279a157612eSMichael Lange DM dmOverlap; 1280a157612eSMichael Lange IS cellPart, part; 1281a157612eSMichael Lange PetscSection cellPartSection, partSection; 128243331d4aSMichael Lange PetscSFNode *remoteRanks, *newRemote; 128343331d4aSMichael Lange const PetscSFNode *oldRemote; 128443331d4aSMichael Lange PetscSF partSF, pointSF, overlapPointSF, overlapSF; 128570034214SMatthew G. Knepley ISLocalToGlobalMapping renumbering; 128670034214SMatthew G. Knepley PetscBool flg; 128770034214SMatthew G. Knepley PetscMPIInt rank, numProcs, p; 128870034214SMatthew G. Knepley PetscErrorCode ierr; 128970034214SMatthew G. Knepley 129070034214SMatthew G. Knepley PetscFunctionBegin; 129170034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 129270034214SMatthew G. Knepley if (sf) PetscValidPointer(sf,4); 129370034214SMatthew G. Knepley PetscValidPointer(dmParallel,5); 129470034214SMatthew G. Knepley 129570034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 129670034214SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 129770034214SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 129870034214SMatthew G. Knepley ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 129970034214SMatthew G. Knepley 130070034214SMatthew G. Knepley *dmParallel = NULL; 130170034214SMatthew G. Knepley if (numProcs == 1) PetscFunctionReturn(0); 130270034214SMatthew G. Knepley 130370034214SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 130470034214SMatthew G. Knepley /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */ 130577623264SMatthew G. Knepley ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 130677623264SMatthew G. Knepley if (overlap > 1) SETERRQ(comm, PETSC_ERR_SUP, "Overlap > 1 not yet implemented"); 130777623264SMatthew G. Knepley ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr); 1308a157612eSMichael Lange ierr = PetscPartitionerPartition(mesh->partitioner, dm, PETSC_FALSE, cellPartSection, &cellPart, NULL, NULL);CHKERRQ(ierr); 130970034214SMatthew G. Knepley /* Create SF assuming a serial partition for all processes: Could check for IS length here */ 131070034214SMatthew G. Knepley if (!rank) numRemoteRanks = numProcs; 131170034214SMatthew G. Knepley else numRemoteRanks = 0; 131270034214SMatthew G. Knepley ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr); 131370034214SMatthew G. Knepley for (p = 0; p < numRemoteRanks; ++p) { 131470034214SMatthew G. Knepley remoteRanks[p].rank = p; 131570034214SMatthew G. Knepley remoteRanks[p].index = 0; 131670034214SMatthew G. Knepley } 131770034214SMatthew G. Knepley ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr); 131870034214SMatthew G. Knepley ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr); 131970034214SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 132070034214SMatthew G. Knepley if (flg) { 1321a157612eSMichael Lange ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr); 132270034214SMatthew G. Knepley ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 132370034214SMatthew G. Knepley ierr = ISView(cellPart, NULL);CHKERRQ(ierr); 132470034214SMatthew G. Knepley ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr); 132570034214SMatthew G. Knepley } 132670034214SMatthew G. Knepley /* Close the partition over the mesh */ 132770034214SMatthew G. Knepley ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr); 132870034214SMatthew G. Knepley /* Create new mesh */ 132970034214SMatthew G. Knepley ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr); 133070034214SMatthew G. Knepley ierr = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr); 133170034214SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr); 133270034214SMatthew G. Knepley pmesh = (DM_Plex*) (*dmParallel)->data; 13334eca1733SMichael Lange pmesh->useAnchors = mesh->useAnchors; 13344eca1733SMichael Lange 133570034214SMatthew G. Knepley /* Distribute sieve points and the global point numbering (replaces creating remote bases) */ 133670034214SMatthew G. Knepley ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr); 133770034214SMatthew G. Knepley if (flg) { 133870034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr); 133970034214SMatthew G. Knepley ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 134070034214SMatthew G. Knepley ierr = ISView(part, NULL);CHKERRQ(ierr); 134170034214SMatthew G. Knepley ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr); 134270034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr); 134370034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr); 134470034214SMatthew G. Knepley } 134577623264SMatthew G. Knepley ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 134670034214SMatthew G. Knepley 1347a157612eSMichael Lange /* Migrate data to a non-overlapping parallel DM */ 1348cc71bff1SMichael Lange ierr = DMPlexDistributeCones(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr); 13490df0e737SMichael Lange ierr = DMPlexDistributeCoordinates(dm, pointSF, *dmParallel);CHKERRQ(ierr); 13502eb1fa14SMichael Lange ierr = DMPlexDistributeLabels(dm, pointSF, *dmParallel);CHKERRQ(ierr); 13519734c634SMichael Lange ierr = DMPlexDistributeSetupHybrid(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr); 1352a6f36705SMichael Lange ierr = DMPlexDistributeSetupTree(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr); 135370034214SMatthew G. Knepley 1354a157612eSMichael Lange /* Build the point SF without overlap */ 1355a157612eSMichael Lange ierr = DMPlexDistributeSF(dm, pointSF, partSection, part, NULL, NULL, *dmParallel);CHKERRQ(ierr); 135670034214SMatthew G. Knepley if (flg) { 135715fff7beSMatthew G. Knepley ierr = PetscPrintf(comm, "Point SF:\n");CHKERRQ(ierr); 135815fff7beSMatthew G. Knepley ierr = PetscSFView((*dmParallel)->sf, NULL);CHKERRQ(ierr); 135970034214SMatthew G. Knepley } 136070034214SMatthew G. Knepley 1361a157612eSMichael Lange if (overlap > 0) { 13623d822a50SMichael Lange ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr); 1363a157612eSMichael Lange /* Add the partition overlap to the distributed DM */ 1364a157612eSMichael Lange ierr = DMPlexDistributeOverlap(*dmParallel, overlap, renumbering, &overlapSF, &dmOverlap);CHKERRQ(ierr); 1365a157612eSMichael Lange ierr = DMDestroy(dmParallel);CHKERRQ(ierr); 1366a157612eSMichael Lange *dmParallel = dmOverlap; 1367c389ea9fSToby Isaac if (flg) { 13683d822a50SMichael Lange ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr); 13693d822a50SMichael Lange ierr = PetscSFView(overlapSF, NULL);CHKERRQ(ierr); 1370c389ea9fSToby Isaac } 137143331d4aSMichael Lange 137243331d4aSMichael Lange /* Re-map the pointSF to establish the full migration pattern */ 137343331d4aSMichael Lange ierr = PetscSFGetGraph(pointSF, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr); 137443331d4aSMichael Lange ierr = PetscSFGetGraph(overlapSF, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr); 137543331d4aSMichael Lange ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr); 137643331d4aSMichael Lange ierr = PetscSFBcastBegin(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr); 137743331d4aSMichael Lange ierr = PetscSFBcastEnd(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr); 137843331d4aSMichael Lange ierr = PetscSFCreate(comm, &overlapPointSF);CHKERRQ(ierr); 137943331d4aSMichael Lange ierr = PetscSFSetGraph(overlapPointSF, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr); 138015fff7beSMatthew G. Knepley ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr); 138143331d4aSMichael Lange ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr); 138243331d4aSMichael Lange pointSF = overlapPointSF; 13833d822a50SMichael Lange ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr); 1384c389ea9fSToby Isaac } 1385bf5244c7SMatthew G. Knepley /* Cleanup Partition */ 1386bf5244c7SMatthew G. Knepley ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr); 1387bf5244c7SMatthew G. Knepley ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr); 1388bf5244c7SMatthew G. Knepley ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr); 1389bf5244c7SMatthew G. Knepley ierr = ISDestroy(&part);CHKERRQ(ierr); 13904eca1733SMichael Lange ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); 13914eca1733SMichael Lange ierr = ISDestroy(&cellPart);CHKERRQ(ierr); 139286719b60SMatthew G. Knepley /* Copy BC */ 139386719b60SMatthew G. Knepley ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr); 139470034214SMatthew G. Knepley /* Cleanup */ 139570034214SMatthew G. Knepley if (sf) {*sf = pointSF;} 139670034214SMatthew G. Knepley else {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);} 139770034214SMatthew G. Knepley ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr); 139870034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 139970034214SMatthew G. Knepley PetscFunctionReturn(0); 140070034214SMatthew G. Knepley } 1401a157612eSMichael Lange 1402a157612eSMichael Lange #undef __FUNCT__ 1403a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap" 1404a157612eSMichael Lange /*@C 1405a157612eSMichael Lange DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM. 1406a157612eSMichael Lange 1407a157612eSMichael Lange Not Collective 1408a157612eSMichael Lange 1409a157612eSMichael Lange Input Parameter: 1410a157612eSMichael Lange + dm - The non-overlapping distrbuted DMPlex object 1411a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default 1412a157612eSMichael Lange 1413a157612eSMichael Lange Output Parameter: 1414a157612eSMichael Lange + sf - The PetscSF used for point distribution 1415a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL 1416a157612eSMichael Lange 1417a157612eSMichael Lange Note: If the mesh was not distributed, the return value is NULL. 1418a157612eSMichael Lange 1419a157612eSMichael Lange The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 1420a157612eSMichael Lange DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 1421a157612eSMichael Lange representation on the mesh. 1422a157612eSMichael Lange 1423a157612eSMichael Lange Level: intermediate 1424a157612eSMichael Lange 1425a157612eSMichael Lange .keywords: mesh, elements 1426a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 1427a157612eSMichael Lange @*/ 1428a157612eSMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, ISLocalToGlobalMapping renumbering, PetscSF *sf, DM *dmOverlap) 1429a157612eSMichael Lange { 1430a157612eSMichael Lange MPI_Comm comm; 1431a157612eSMichael Lange PetscMPIInt rank; 1432a157612eSMichael Lange PetscSection rootSection, leafSection; 1433a157612eSMichael Lange IS rootrank, leafrank; 1434a157612eSMichael Lange PetscSection coneSection; 1435a157612eSMichael Lange PetscSF overlapSF, migrationSF, pointSF, newPointSF; 1436a157612eSMichael Lange PetscSFNode *ghostRemote; 1437a157612eSMichael Lange const PetscSFNode *overlapRemote; 1438a157612eSMichael Lange ISLocalToGlobalMapping overlapRenumbering; 1439a157612eSMichael Lange const PetscInt *renumberingArray, *overlapLocal; 1440a157612eSMichael Lange PetscInt dim, p, pStart, pEnd, conesSize, idx; 1441a157612eSMichael Lange PetscInt numGhostPoints, numOverlapPoints, numSharedPoints, overlapLeaves; 1442a157612eSMichael Lange PetscInt *cones, *ghostLocal, *overlapRenumberingArray, *pointIDs, *recvPointIDs; 1443a157612eSMichael Lange PetscErrorCode ierr; 1444a157612eSMichael Lange 1445a157612eSMichael Lange PetscFunctionBegin; 1446a157612eSMichael Lange PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1447a157612eSMichael Lange if (sf) PetscValidPointer(sf, 3); 1448a157612eSMichael Lange PetscValidPointer(dmOverlap, 4); 1449a157612eSMichael Lange 1450a157612eSMichael Lange ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1451a157612eSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1452a157612eSMichael Lange ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1453a157612eSMichael Lange 1454a157612eSMichael Lange /* Compute point overlap with neighbouring processes on the distributed DM */ 1455a157612eSMichael Lange ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1456a157612eSMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 1457a157612eSMichael Lange ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); 1458a157612eSMichael Lange ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr); 1459a157612eSMichael Lange ierr = DMPlexCreateOverlap(dm, rootSection, rootrank, leafSection, leafrank, &overlapSF);CHKERRQ(ierr); 146015fff7beSMatthew G. Knepley ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 146115fff7beSMatthew G. Knepley ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 146215fff7beSMatthew G. Knepley ierr = ISDestroy(&rootrank);CHKERRQ(ierr); 146315fff7beSMatthew G. Knepley ierr = ISDestroy(&leafrank);CHKERRQ(ierr); 1464a157612eSMichael Lange ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1465a157612eSMichael Lange 1466a157612eSMichael Lange /* Build dense migration SF that maps the non-overlapping partition to the overlapping one */ 1467a157612eSMichael Lange ierr = DMPlexCreateOverlapMigrationSF(dm, overlapSF, &migrationSF);CHKERRQ(ierr); 1468a157612eSMichael Lange 1469a157612eSMichael Lange /* Convert cones to global numbering before migrating them */ 1470a157612eSMichael Lange ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr); 1471a157612eSMichael Lange ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr); 1472a157612eSMichael Lange ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 1473a157612eSMichael Lange ierr = ISLocalToGlobalMappingApplyBlock(renumbering, conesSize, cones, cones);CHKERRQ(ierr); 1474a157612eSMichael Lange 1475a157612eSMichael Lange /* Derive the new local-to-global mapping from the old one */ 1476a157612eSMichael Lange ierr = PetscSFGetGraph(migrationSF, NULL, &overlapLeaves, &overlapLocal, &overlapRemote);CHKERRQ(ierr); 1477a157612eSMichael Lange ierr = PetscMalloc1(overlapLeaves, &overlapRenumberingArray);CHKERRQ(ierr); 1478a157612eSMichael Lange ierr = ISLocalToGlobalMappingGetBlockIndices(renumbering, &renumberingArray);CHKERRQ(ierr); 1479729b3788SMatthew G. Knepley ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, (PetscInt *) renumberingArray, overlapRenumberingArray);CHKERRQ(ierr); 1480729b3788SMatthew G. Knepley ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, (PetscInt *) renumberingArray, overlapRenumberingArray);CHKERRQ(ierr); 1481a157612eSMichael Lange ierr = ISLocalToGlobalMappingCreate(comm, 1, overlapLeaves, (const PetscInt*) overlapRenumberingArray, PETSC_OWN_POINTER, &overlapRenumbering);CHKERRQ(ierr); 1482a157612eSMichael Lange 1483a157612eSMichael Lange /* Build the overlapping DM */ 1484a157612eSMichael Lange ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr); 1485a157612eSMichael Lange ierr = DMSetDimension(*dmOverlap, dim);CHKERRQ(ierr); 1486a157612eSMichael Lange ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr); 1487a157612eSMichael Lange ierr = DMPlexDistributeCones(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr); 1488a157612eSMichael Lange ierr = DMPlexDistributeCoordinates(dm, migrationSF, *dmOverlap);CHKERRQ(ierr); 1489a157612eSMichael Lange ierr = DMPlexDistributeLabels(dm, migrationSF, *dmOverlap);CHKERRQ(ierr); 1490a157612eSMichael Lange ierr = DMPlexDistributeSetupHybrid(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr); 1491a157612eSMichael Lange 1492a157612eSMichael Lange /* Build the new point SF by propagating the depthShift generate remote root indices */ 1493a157612eSMichael Lange ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr); 1494a157612eSMichael Lange ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, NULL, NULL);CHKERRQ(ierr); 1495a157612eSMichael Lange ierr = PetscSFGetGraph(overlapSF, NULL, &numOverlapPoints, NULL, NULL);CHKERRQ(ierr); 1496a157612eSMichael Lange numGhostPoints = numSharedPoints + numOverlapPoints; 149709b7985cSMatthew G. Knepley ierr = PetscMalloc1(numGhostPoints, &ghostLocal);CHKERRQ(ierr); 149809b7985cSMatthew G. Knepley ierr = PetscMalloc1(numGhostPoints, &ghostRemote);CHKERRQ(ierr); 1499a157612eSMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1500a157612eSMichael Lange ierr = PetscMalloc2(pEnd-pStart, &pointIDs, overlapLeaves, &recvPointIDs);CHKERRQ(ierr); 1501a157612eSMichael Lange for (p=0; p<overlapLeaves; p++) { 1502a157612eSMichael Lange if (overlapRemote[p].rank == rank) pointIDs[overlapRemote[p].index] = overlapLocal[p]; 1503a157612eSMichael Lange } 1504a157612eSMichael Lange ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr); 1505a157612eSMichael Lange ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr); 1506a157612eSMichael Lange for (idx=0, p=0; p<overlapLeaves; p++) { 1507a157612eSMichael Lange if (overlapRemote[p].rank != rank) { 1508a157612eSMichael Lange ghostLocal[idx] = overlapLocal[p]; 1509a157612eSMichael Lange ghostRemote[idx].index = recvPointIDs[p]; 1510a157612eSMichael Lange ghostRemote[idx].rank = overlapRemote[p].rank; 1511a157612eSMichael Lange idx++; 1512a157612eSMichael Lange } 1513a157612eSMichael Lange } 1514a157612eSMichael Lange ierr = DMPlexGetChart(*dmOverlap, &pStart, &pEnd);CHKERRQ(ierr); 1515a157612eSMichael Lange ierr = PetscSFCreate(comm, &newPointSF);;CHKERRQ(ierr); 1516a157612eSMichael Lange ierr = PetscSFSetGraph(newPointSF, pEnd - pStart, numGhostPoints, ghostLocal, PETSC_OWN_POINTER, ghostRemote, PETSC_OWN_POINTER);CHKERRQ(ierr); 1517a157612eSMichael Lange ierr = DMSetPointSF(*dmOverlap, newPointSF);CHKERRQ(ierr); 151815fff7beSMatthew G. Knepley ierr = PetscSFDestroy(&newPointSF);CHKERRQ(ierr); 1519a157612eSMichael Lange /* Cleanup overlap partition */ 1520a157612eSMichael Lange ierr = ISLocalToGlobalMappingDestroy(&overlapRenumbering);CHKERRQ(ierr); 1521a157612eSMichael Lange ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr); 1522a157612eSMichael Lange ierr = PetscFree2(pointIDs, recvPointIDs);CHKERRQ(ierr); 1523a157612eSMichael Lange if (sf) *sf = migrationSF; 1524a157612eSMichael Lange else {ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr);} 1525a157612eSMichael Lange PetscFunctionReturn(0); 1526a157612eSMichael Lange } 1527