/*@ DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution Collective on DM Input Parameters: + dm - The DMPlex object . pointSF - The PetscSF describing the communication pattern . originalSection - The PetscSection for existing data layout - originalVec - The existing data Output Parameters: + newSection - The PetscSF describing the new data layout - newVec - The new data Level: developer .seealso: DMPlexDistribute(), DMPlexDistributeData() @*/ PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec) { PetscSF fieldSF; PetscInt *remoteOffsets, fieldSize; PetscScalar *originalValues, *newValues; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr); ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr); ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr); ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr); ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr); ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode DistributeCoordinates(DM dm, PetscSF pointSF, DM parallelDM) { PetscSF coordSF; PetscSection originalCoordSection, newCoordSection; Vec coordinates, newCoordinates; PetscScalar *coords, *newCoords; PetscInt *remoteOffsets, coordSize; PetscErrorCode ierr; PetscFunctionBegin; ierr = DMMeshGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); ierr = DMMeshGetCoordinateSection(parallelDM, &newCoordSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(pointSF, originalCoordSection, &remoteOffsets, newCoordSection);CHKERRQ(ierr); ierr = DMMeshGetCoordinateVec(dm, &coordinates);CHKERRQ(ierr); ierr = DMMeshGetCoordinateVec(parallelDM, &newCoordinates);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(newCoordSection, &coordSize);CHKERRQ(ierr); ierr = VecSetSizes(newCoordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); ierr = VecSetFromOptions(newCoordinates);CHKERRQ(ierr); ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); ierr = VecGetArray(newCoordinates, &newCoords);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(pointSF, originalCoordSection, remoteOffsets, newCoordSection, &coordSF);CHKERRQ(ierr); ierr = PetscSFBcastBegin(coordSF, MPIU_SCALAR, coords, newCoords);CHKERRQ(ierr); ierr = PetscSFBcastEnd(coordSF, MPIU_SCALAR, coords, newCoords);CHKERRQ(ierr); ierr = PetscSFDestroy(&coordSF);CHKERRQ(ierr); ierr = VecRestoreArray(newCoordinates, &newCoords);CHKERRQ(ierr); ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode DMLabelDistribute_Internal(DMLabel label, PetscSF sf, PetscSection *leafSection, PetscInt **leafStrata) { MPI_Comm comm; PetscInt s, l, nroots, nleaves, dof, offset, size; PetscInt *remoteOffsets, *rootStrata, *rootIdx; PetscSection rootSection; PetscSF labelSF; PetscErrorCode ierr; PetscFunctionBegin; if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); /* Build a section of stratum values per point, generate the according SF and distribute point-wise stratum values to leaves. */ ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); if (label) { for (s = 0; s < label->numStrata; ++s) { for (l = 0; l < label->stratumSizes[s]; l++) { ierr = PetscSectionGetDof(rootSection, label->points[s][l], &dof);CHKERRQ(ierr); ierr = PetscSectionSetDof(rootSection, label->points[s][l], dof+1);CHKERRQ(ierr); } } } ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); /* Create a point-wise array of stratum values */ ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); if (label) { for (s = 0; s < label->numStrata; ++s) { for (l = 0; l < label->stratumSizes[s]; l++) { const PetscInt p = label->points[s][l]; ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); rootStrata[offset+rootIdx[p]++] = label->stratumValues[s]; } } } /* Build SF that maps label points to remote processes */ ierr = PetscSectionCreate(comm, leafSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, *leafSection);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, *leafSection, &labelSF);CHKERRQ(ierr); ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); /* Send the strata for each point over the derived SF */ ierr = PetscSectionGetStorageSize(*leafSection, &size);CHKERRQ(ierr); ierr = PetscMalloc1(size, leafStrata);CHKERRQ(ierr); ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, *leafStrata);CHKERRQ(ierr); /* Clean up */ ierr = PetscFree(rootStrata);CHKERRQ(ierr); ierr = PetscFree(rootIdx);CHKERRQ(ierr); ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution Collective on DM Input Parameters: + dm - The DMPlex object . pointSF - The PetscSF describing the communication pattern . originalSection - The PetscSection for existing data layout . datatype - The type of data - originalData - The existing data Output Parameters: + newSection - The PetscSF describing the new data layout - newData - The new data Level: developer .seealso: DMPlexDistribute(), DMPlexDistributeField() @*/ PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData) { PetscSF fieldSF; PetscInt *remoteOffsets, fieldSize; PetscMPIInt dataSize; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr); ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ DMNetworkDistribute - Distributes the network and moves associated component data. Collective Input Parameter: + oldDM - the original DMNetwork object - overlap - The overlap of partitions, 0 is the default Output Parameter: . distDM - the distributed DMNetwork object Notes: This routine should be called only when using multiple processors. Distributes the network with <overlap>-overlapping partitioning of the edges. Level: intermediate .seealso: DMNetworkCreate @*/ PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM) { PetscErrorCode ierr; DM_Network *oldDMnetwork = (DM_Network*)oldDM->data; PetscSF pointsf; DM newDM; DM_Network *newDMnetwork; PetscFunctionBegin; ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr); newDMnetwork = (DM_Network*)newDM->data; newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); /* Distribute plex dm and dof section */ ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); /* Distribute dof section */ ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr); /* Distribute data and associated section */ ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPI_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); /* Destroy point SF */ ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart; newDMnetwork->NNodes = oldDMnetwork->NNodes; newDMnetwork->NEdges = oldDMnetwork->NEdges; /* Set Dof section as the default section for dm */ ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); *distDM = newDM; PetscFunctionReturn(0); }
/*@C DMPlexDistribute - Distributes the mesh and any associated sections. Not Collective Input Parameter: + dm - The original DMPlex object . partitioner - The partitioning package, or NULL for the default - overlap - The overlap of partitions, 0 is the default Output Parameter: + sf - The PetscSF used for point distribution - parallelMesh - The distributed DMPlex object, or NULL Note: If the mesh was not distributed, the return value is NULL. The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function representation on the mesh. Level: intermediate .keywords: mesh, elements .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() @*/ PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel) { DM_Plex *mesh = (DM_Plex*) dm->data, *pmesh; MPI_Comm comm; const PetscInt height = 0; PetscInt dim, numRemoteRanks; IS origCellPart, origPart, cellPart, part; PetscSection origCellPartSection, origPartSection, cellPartSection, partSection; PetscSFNode *remoteRanks; PetscSF partSF, pointSF, coneSF; ISLocalToGlobalMapping renumbering; PetscSection originalConeSection, newConeSection; PetscInt *remoteOffsets; PetscInt *cones, *newCones, newConesSize; PetscBool flg; PetscMPIInt rank, numProcs, p; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(dm, DM_CLASSID, 1); if (sf) PetscValidPointer(sf,4); PetscValidPointer(dmParallel,5); ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); *dmParallel = NULL; if (numProcs == 1) PetscFunctionReturn(0); ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */ ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr); if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented"); ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr); /* Create SF assuming a serial partition for all processes: Could check for IS length here */ if (!rank) numRemoteRanks = numProcs; else numRemoteRanks = 0; ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr); for (p = 0; p < numRemoteRanks; ++p) { remoteRanks[p].rank = p; remoteRanks[p].index = 0; } ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr); ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); if (flg) { ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr); ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISView(cellPart, NULL);CHKERRQ(ierr); if (origCellPart) { ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr); ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISView(origCellPart, NULL);CHKERRQ(ierr); } ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr); } /* Close the partition over the mesh */ ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr); ierr = ISDestroy(&cellPart);CHKERRQ(ierr); ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); /* Create new mesh */ ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr); ierr = DMPlexSetDimension(*dmParallel, dim);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr); pmesh = (DM_Plex*) (*dmParallel)->data; /* Distribute sieve points and the global point numbering (replaces creating remote bases) */ ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr); if (flg) { ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr); ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISView(part, NULL);CHKERRQ(ierr); ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr); ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr); ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr); } ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr); ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); /* Distribute cone section */ ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); ierr = DMSetUp(*dmParallel);CHKERRQ(ierr); { PetscInt pStart, pEnd, p; ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr); for (p = pStart; p < pEnd; ++p) { PetscInt coneSize; ierr = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr); pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize); } } /* Communicate and renumber cones */ ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr); ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr); ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr); if (flg) { ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr); } ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr); ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr); ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); /* Create supports and stratify sieve */ { PetscInt pStart, pEnd; ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr); ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr); } ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr); ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr); /* Distribute Coordinates */ { PetscSection originalCoordSection, newCoordSection; Vec originalCoordinates, newCoordinates; PetscInt bs; const char *name; ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr); ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr); ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr); ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr); ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr); ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr); ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr); ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr); ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); } /* Distribute labels */ ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); { DMLabel next = mesh->labels, newNext = pmesh->labels; PetscInt numLabels = 0, l; /* Bcast number of labels */ while (next) {++numLabels; next = next->next;} ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); next = mesh->labels; for (l = 0; l < numLabels; ++l) { DMLabel labelNew; PetscBool isdepth; /* Skip "depth" because it is recreated */ if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);} ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); if (isdepth) {if (!rank) next = next->next; continue;} ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr); /* Insert into list */ if (newNext) newNext->next = labelNew; else pmesh->labels = labelNew; newNext = labelNew; if (!rank) next = next->next; } } ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); /* Setup hybrid structure */ { const PetscInt *gpoints; PetscInt depth, n, d; for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];} ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr); ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); for (d = 0; d <= dim; ++d) { PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p; if (pmax < 0) continue; ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr); ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr); ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr); for (p = 0; p < n; ++p) { const PetscInt point = gpoints[p]; if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax; } if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax; else pmesh->hybridPointMax[d] = -1; } ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr); } /* Cleanup Partition */ ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr); ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr); ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr); ierr = ISDestroy(&part);CHKERRQ(ierr); /* Create point SF for parallel mesh */ ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); { const PetscInt *leaves; PetscSFNode *remotePoints, *rowners, *lowners; PetscInt numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints; PetscInt pStart, pEnd; ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr); ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr); ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr); for (p=0; p<numRoots; p++) { rowners[p].rank = -1; rowners[p].index = -1; } if (origCellPart) { /* Make sure points in the original partition are not assigned to other procs */ const PetscInt *origPoints; ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr); ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr); for (p = 0; p < numProcs; ++p) { PetscInt dof, off, d; ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr); ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr); for (d = off; d < off+dof; ++d) { rowners[origPoints[d]].rank = p; } } ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr); ierr = ISDestroy(&origPart);CHKERRQ(ierr); ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr); } ierr = ISDestroy(&origCellPart);CHKERRQ(ierr); ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr); ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); for (p = 0; p < numLeaves; ++p) { if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */ lowners[p].rank = rank; lowners[p].index = leaves ? leaves[p] : p; } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */ lowners[p].rank = -2; lowners[p].index = -2; } } for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */ rowners[p].rank = -3; rowners[p].index = -3; } ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); for (p = 0; p < numLeaves; ++p) { if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed"); if (lowners[p].rank != rank) ++numGhostPoints; } ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr); ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr); for (p = 0, gp = 0; p < numLeaves; ++p) { if (lowners[p].rank != rank) { ghostPoints[gp] = leaves ? leaves[p] : p; remotePoints[gp].rank = lowners[p].rank; remotePoints[gp].index = lowners[p].index; ++gp; } } ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr); ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr); } pmesh->useCone = mesh->useCone; pmesh->useClosure = mesh->useClosure; ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); /* Copy BC */ ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr); /* Cleanup */ if (sf) {*sf = pointSF;} else {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);} ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr); ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) { MPI_Comm comm; PetscSection rootSection, leafSection; PetscSF labelSF; PetscInt p, pStart, pEnd, l, lStart, lEnd, s, nroots, nleaves, size, dof, offset, stratum; PetscInt *remoteOffsets, *rootStrata, *rootIdx, *leafStrata, *strataIdx; char *name; PetscInt nameSize; size_t len = 0; PetscMPIInt rank, numProcs; PetscErrorCode ierr; PetscFunctionBegin; if (label) {ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr);} ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); /* Bcast name */ if (!rank) {ierr = PetscStrlen(label->name, &len);CHKERRQ(ierr);} nameSize = len; ierr = MPI_Bcast(&nameSize, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); ierr = PetscMalloc1(nameSize+1, &name);CHKERRQ(ierr); if (!rank) {ierr = PetscMemcpy(name, label->name, nameSize+1);CHKERRQ(ierr);} ierr = MPI_Bcast(name, nameSize+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); ierr = DMLabelCreate(name, labelNew);CHKERRQ(ierr); ierr = PetscFree(name);CHKERRQ(ierr); /* Bcast numStrata */ if (!rank) (*labelNew)->numStrata = label->numStrata; ierr = MPI_Bcast(&(*labelNew)->numStrata, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); /* Bcast stratumValues */ ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); if (!rank) {ierr = PetscMemcpy((*labelNew)->stratumValues, label->stratumValues, label->numStrata * sizeof(PetscInt));CHKERRQ(ierr);} ierr = MPI_Bcast((*labelNew)->stratumValues, (*labelNew)->numStrata, MPIU_INT, 0, comm);CHKERRQ(ierr); ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->arrayValid);CHKERRQ(ierr); for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->arrayValid[s] = PETSC_TRUE; /* Build a section detailing strata-per-point, distribute and build SF from that and then send our points. */ ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); ierr = PetscSectionSetChart(rootSection, 0, nroots);CHKERRQ(ierr); if (label) { for (s = 0; s < label->numStrata; ++s) { lStart = 0; lEnd = label->stratumSizes[s]; for (l=lStart; l<lEnd; l++) { ierr = PetscSectionGetDof(rootSection, label->points[s][l], &dof);CHKERRQ(ierr); ierr = PetscSectionSetDof(rootSection, label->points[s][l], dof+1);CHKERRQ(ierr); } } } ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); /* Create a point-wise array of point strata */ ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); ierr = PetscMalloc1(size, &rootStrata);CHKERRQ(ierr); ierr = PetscCalloc1(nroots, &rootIdx);CHKERRQ(ierr); if (label) { for (s = 0; s < label->numStrata; ++s) { lStart = 0; lEnd = label->stratumSizes[s]; for (l=lStart; l<lEnd; l++) { p = label->points[s][l]; ierr = PetscSectionGetOffset(rootSection, p, &offset);CHKERRQ(ierr); rootStrata[offset+rootIdx[p]++] = s; } } } /* Build SF that maps label points to remote processes */ ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(sf, rootSection, &remoteOffsets, leafSection);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(sf, rootSection, remoteOffsets, leafSection, &labelSF);CHKERRQ(ierr); /* Send the strata for each point over the derived SF */ ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); ierr = PetscMalloc1(size, &leafStrata);CHKERRQ(ierr); ierr = PetscSFBcastBegin(labelSF, MPIU_INT, rootStrata, leafStrata);CHKERRQ(ierr); ierr = PetscSFBcastEnd(labelSF, MPIU_INT, rootStrata, leafStrata);CHKERRQ(ierr); /* Rebuild the point strata on the receiver */ ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->stratumSizes);CHKERRQ(ierr); ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); for (p=pStart; p<pEnd; p++) { ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); for (s=0; s<dof; s++) { (*labelNew)->stratumSizes[leafStrata[offset+s]]++; } } ierr = PetscCalloc1((*labelNew)->numStrata,&(*labelNew)->ht);CHKERRQ(ierr); ierr = PetscMalloc1((*labelNew)->numStrata,&(*labelNew)->points);CHKERRQ(ierr); for (s = 0; s < (*labelNew)->numStrata; ++s) { PetscHashICreate((*labelNew)->ht[s]); ierr = PetscMalloc1((*labelNew)->stratumSizes[s], &(*labelNew)->points[s]);CHKERRQ(ierr); } /* Insert points into new strata */ ierr = PetscCalloc1((*labelNew)->numStrata, &strataIdx);CHKERRQ(ierr); ierr = PetscSectionGetChart(leafSection, &pStart, &pEnd);CHKERRQ(ierr); for (p=pStart; p<pEnd; p++) { ierr = PetscSectionGetDof(leafSection, p, &dof);CHKERRQ(ierr); ierr = PetscSectionGetOffset(leafSection, p, &offset);CHKERRQ(ierr); for (s=0; s<dof; s++) { stratum = leafStrata[offset+s]; (*labelNew)->points[stratum][strataIdx[stratum]++] = p; } } ierr = PetscFree(rootStrata);CHKERRQ(ierr); ierr = PetscFree(leafStrata);CHKERRQ(ierr); ierr = PetscFree(rootIdx);CHKERRQ(ierr); ierr = PetscFree(strataIdx);CHKERRQ(ierr); ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); ierr = PetscSFDestroy(&labelSF);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* Distribute cones - Partitioning: input partition point map and naive sf, output sf with inverse of map, distribute points - Distribute section: input current sf, communicate sizes and offsets, output local section and offsets (only use for new sf) - Create SF for values: input current sf and offsets, output new sf - Distribute values: input new sf, communicate values */ PetscErrorCode DistributeMesh(DM dm, AppCtx *user, PetscSF *pointSF, DM *parallelDM) { MPI_Comm comm = ((PetscObject) dm)->comm; const PetscInt height = 0; PetscInt dim, numRemoteRanks; IS cellPart, part; PetscSection cellPartSection, partSection; PetscSFNode *remoteRanks; PetscSF partSF; ISLocalToGlobalMapping renumbering; PetscSF coneSF; PetscSection originalConeSection, newConeSection; PetscInt *remoteOffsets, newConesSize; PetscInt *cones, *newCones; PetscMPIInt numProcs, rank, p; PetscErrorCode ierr; PetscFunctionBegin; ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = DMMeshGetDimension(dm, &dim);CHKERRQ(ierr); /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */ ierr = DMMeshCreatePartition(dm, &cellPartSection, &cellPart, height);CHKERRQ(ierr); /* Create SF assuming a serial partition for all processes: Could check for IS length here */ if (!rank) { numRemoteRanks = numProcs; } else { numRemoteRanks = 0; } ierr = PetscMalloc(numRemoteRanks * sizeof(PetscSFNode), &remoteRanks);CHKERRQ(ierr); for(p = 0; p < numRemoteRanks; ++p) { remoteRanks[p].rank = p; remoteRanks[p].index = 0; } ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr); ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, PETSC_NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr); /* Debugging */ ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr); ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISView(cellPart, PETSC_NULL);CHKERRQ(ierr); ierr = PetscSFView(partSF, PETSC_NULL);CHKERRQ(ierr); /* Close the partition over the mesh */ ierr = DMMeshCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr); ierr = ISDestroy(&cellPart);CHKERRQ(ierr); ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); /* Create new mesh */ ierr = DMMeshCreate(comm, parallelDM);CHKERRQ(ierr); ierr = DMMeshSetDimension(*parallelDM, dim);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) *parallelDM, "Parallel Mesh");CHKERRQ(ierr); /* Distribute sieve points and the global point numbering (replaces creating remote bases) */ ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, pointSF);CHKERRQ(ierr); /* Debugging */ ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr); ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = ISView(part, PETSC_NULL);CHKERRQ(ierr); ierr = PetscSFView(*pointSF, PETSC_NULL);CHKERRQ(ierr); ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr); ierr = ISLocalToGlobalMappingView(renumbering, PETSC_NULL);CHKERRQ(ierr); /* Cleanup */ ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr); ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr); ierr = ISDestroy(&part);CHKERRQ(ierr); /* Distribute cone section */ ierr = DMMeshGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); ierr = DMMeshGetConeSection(*parallelDM, &newConeSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(*pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); ierr = DMMeshSetUp(*parallelDM);CHKERRQ(ierr); /* Communicate and renumber cones */ ierr = PetscSFCreateSectionSF(*pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); ierr = DMMeshGetCones(dm, &cones);CHKERRQ(ierr); ierr = DMMeshGetCones(*parallelDM, &newCones);CHKERRQ(ierr); ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); ierr = ISGlobalToLocalMappingApply(renumbering, IS_GTOLM_MASK, newConesSize, newCones, PETSC_NULL, newCones);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr); /* Debugging */ ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscSFView(coneSF, PETSC_NULL);CHKERRQ(ierr); ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); /* Create supports and stratify sieve */ ierr = DMMeshSymmetrize(*parallelDM);CHKERRQ(ierr); ierr = DMMeshStratify(*parallelDM);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* This interpolates the PointSF in parallel following local interpolation */ static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth) { PetscMPIInt numProcs, rank; PetscInt p, c, d, dof, offset; PetscInt numLeaves, numRoots, candidatesSize, candidatesRemoteSize; const PetscInt *localPoints; const PetscSFNode *remotePoints; PetscSFNode *candidates, *candidatesRemote, *claims; PetscSection candidateSection, candidateSectionRemote, claimSection; PetscHashI leafhash; PetscHashIJ roothash; PetscHashIJKey key; PetscErrorCode ierr; PetscFunctionBegin; ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); if (numProcs < 2 || numRoots < 0) PetscFunctionReturn(0); /* Build hashes of points in the SF for efficient lookup */ PetscHashICreate(leafhash); PetscHashIJCreate(&roothash); ierr = PetscHashIJSetMultivalued(roothash, PETSC_FALSE);CHKERRQ(ierr); for (p = 0; p < numLeaves; ++p) { PetscHashIAdd(leafhash, localPoints[p], p); key.i = remotePoints[p].index; key.j = remotePoints[p].rank; PetscHashIJAdd(roothash, key, p); } /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves, where each candidate is defined by the root entry for the other vertex that defines the edge. */ ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); { PetscInt leaf, root, idx, a, *adj = NULL; for (p = 0; p < numLeaves; ++p) { PetscInt adjSize = PETSC_DETERMINE; ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); for (a = 0; a < adjSize; ++a) { PetscHashIMap(leafhash, adj[a], leaf); if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);} } } ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); for (p = 0; p < numLeaves; ++p) { PetscInt adjSize = PETSC_DETERMINE; ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr); ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); for (idx = 0, a = 0; a < adjSize; ++a) { PetscHashIMap(leafhash, adj[a], root); if (root >= 0) candidates[offset+idx++] = remotePoints[root]; } } ierr = PetscFree(adj);CHKERRQ(ierr); } /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ { PetscSF sfMulti, sfInverse, sfCandidates; PetscInt *remoteOffsets; ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); } /* Walk local roots and check for each remote candidate whether we know all required points, either from owning it or having a root entry in the point SF. If we do we place a claim by replacing the vertex number with our edge ID. */ { PetscInt idx, root, joinSize, vertices[2]; const PetscInt *rootdegree, *join = NULL; ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); /* Loop remote edge connections and put in a claim if both vertices are known */ for (idx = 0, p = 0; p < numRoots; ++p) { for (d = 0; d < rootdegree[p]; ++d) { ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); for (c = 0; c < dof; ++c) { /* We own both vertices, so we claim the edge by replacing vertex with edge */ if (candidatesRemote[offset+c].rank == rank) { vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index; ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); if (joinSize == 1) candidatesRemote[offset+c].index = join[0]; ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); continue; } /* If we own one vertex and share a root with the other, we claim it */ key.i = candidatesRemote[offset+c].index; key.j = candidatesRemote[offset+c].rank; PetscHashIJGet(roothash, key, &root); if (root >= 0) { vertices[0] = p; vertices[1] = localPoints[root]; ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); if (joinSize == 1) { candidatesRemote[offset+c].index = join[0]; candidatesRemote[offset+c].rank = rank; } ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); } } idx++; } } } /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ { PetscSF sfMulti, sfClaims, sfPointNew; PetscHashI claimshash; PetscInt size, pStart, pEnd, root, joinSize, numLocalNew; PetscInt *remoteOffsets, *localPointsNew, vertices[2]; const PetscInt *join = NULL; PetscSFNode *remotePointsNew; ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr); ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr); ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); /* Walk the original section of local supports and add an SF entry for each updated item */ PetscHashICreate(claimshash); for (p = 0; p < numRoots; ++p) { ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); for (d = 0; d < dof; ++d) { if (candidates[offset+d].index != claims[offset+d].index) { key.i = candidates[offset+d].index; key.j = candidates[offset+d].rank; PetscHashIJGet(roothash, key, &root); if (root >= 0) { vertices[0] = p; vertices[1] = localPoints[root]; ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); if (joinSize == 1) PetscHashIAdd(claimshash, join[0], offset+d); ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); } } } } /* Create new pointSF from hashed claims */ PetscHashISize(claimshash, numLocalNew); ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); for (p = 0; p < numLeaves; ++p) { localPointsNew[p] = localPoints[p]; remotePointsNew[p].index = remotePoints[p].index; remotePointsNew[p].rank = remotePoints[p].rank; } p = numLeaves; ierr = PetscHashIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { PetscHashIMap(claimshash, localPointsNew[p], offset); remotePointsNew[p] = claims[offset]; } ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); PetscHashIDestroy(claimshash); } PetscHashIDestroy(leafhash); ierr = PetscHashIJDestroy(&roothash);CHKERRQ(ierr); ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); ierr = PetscFree(candidates);CHKERRQ(ierr); ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); ierr = PetscFree(claims);CHKERRQ(ierr); PetscFunctionReturn(0); }