/* DMLabelMakeValid_Private - Transfer stratum data from the hash format to the sorted list format Input parameter: + label - The DMLabel - v - The stratum value Output parameter: . label - The DMLabel with stratum in sorted list format Level: developer .seealso: DMLabelCreate() */ static PetscErrorCode DMLabelMakeValid_Private(DMLabel label, PetscInt v) { PetscInt off; PetscErrorCode ierr; if (label->arrayValid[v]) return 0; if (v >= label->numStrata) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Trying to access invalid stratum %D in DMLabelMakeValid_Private\n", v); PetscFunctionBegin; PetscHashISize(label->ht[v], label->stratumSizes[v]); ierr = PetscMalloc1(label->stratumSizes[v], &label->points[v]);CHKERRQ(ierr); off = 0; ierr = PetscHashIGetKeys(label->ht[v], &off, &(label->points[v][0]));CHKERRQ(ierr); if (off != label->stratumSizes[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of contributed points %D from value %D should be %D", off, label->stratumValues[v], label->stratumSizes[v]); PetscHashIClear(label->ht[v]); ierr = PetscSortInt(label->stratumSizes[v], label->points[v]);CHKERRQ(ierr); if (label->bt) { PetscInt p; for (p = 0; p < label->stratumSizes[v]; ++p) { const PetscInt point = label->points[v][p]; if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); ierr = PetscBTSet(label->bt, point - label->pStart);CHKERRQ(ierr); } } label->arrayValid[v] = PETSC_TRUE; ++label->state; PetscFunctionReturn(0); }
PetscErrorCode DMLabelDistribute(DMLabel label, PetscSF sf, DMLabel *labelNew) { MPI_Comm comm; PetscSection leafSection; PetscInt p, pStart, pEnd, s, size, dof, offset, stratum; PetscInt *leafStrata, *strataIdx; char *name; PetscInt nameSize; PetscHashI stratumHash; PETSC_UNUSED PetscHashIIter ret, iter; size_t len = 0; PetscMPIInt rank; 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); /* 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 defaultValue */ if (!rank) (*labelNew)->defaultValue = label->defaultValue; ierr = MPI_Bcast(&(*labelNew)->defaultValue, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); /* Distribute stratum values over the SF and get the point mapping on the receiver */ ierr = DMLabelDistribute_Internal(label, sf, &leafSection, &leafStrata);CHKERRQ(ierr); /* Determine received stratum values and initialise new label*/ PetscHashICreate(stratumHash); ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); for (p = 0; p < size; ++p) PetscHashIPut(stratumHash, leafStrata[p], ret, iter); PetscHashISize(stratumHash, (*labelNew)->numStrata); ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->arrayValid);CHKERRQ(ierr); for (s = 0; s < (*labelNew)->numStrata; ++s) (*labelNew)->arrayValid[s] = PETSC_TRUE; ierr = PetscMalloc1((*labelNew)->numStrata, &(*labelNew)->stratumValues);CHKERRQ(ierr); /* Turn leafStrata into indices rather than stratum values */ offset = 0; ierr = PetscHashIGetKeys(stratumHash, &offset, (*labelNew)->stratumValues);CHKERRQ(ierr); for (s = 0; s < (*labelNew)->numStrata; ++s) { for (p = 0; p < size; ++p) { if (leafStrata[p] == (*labelNew)->stratumValues[s]) leafStrata[p] = s; } } /* 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; } } PetscHashIDestroy(stratumHash); ierr = PetscFree(leafStrata);CHKERRQ(ierr); ierr = PetscFree(strataIdx);CHKERRQ(ierr); ierr = PetscSectionDestroy(&leafSection);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); }