1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/matorderimpl.h> /*I "petscmat.h" I*/ 3 4 static PetscErrorCode DMPlexCreateOrderingClosure_Static(DM dm, PetscInt numPoints, const PetscInt pperm[], PetscInt **clperm, PetscInt **invclperm) 5 { 6 PetscInt *perm, *iperm; 7 PetscInt depth, d, pStart, pEnd, fStart, fMax, fEnd, p; 8 PetscErrorCode ierr; 9 10 PetscFunctionBegin; 11 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 12 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 13 ierr = PetscMalloc1(pEnd-pStart,&perm);CHKERRQ(ierr); 14 ierr = PetscMalloc1(pEnd-pStart,&iperm);CHKERRQ(ierr); 15 for (p = pStart; p < pEnd; ++p) iperm[p] = -1; 16 for (d = depth; d > 0; --d) { 17 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 18 ierr = DMPlexGetDepthStratum(dm, d-1, &fStart, &fEnd);CHKERRQ(ierr); 19 fMax = fStart; 20 for (p = pStart; p < pEnd; ++p) { 21 const PetscInt *cone; 22 PetscInt point, coneSize, c; 23 24 if (d == depth) { 25 perm[p] = pperm[p]; 26 iperm[pperm[p]] = p; 27 } 28 point = perm[p]; 29 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 30 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 31 for (c = 0; c < coneSize; ++c) { 32 const PetscInt oldc = cone[c]; 33 const PetscInt newc = iperm[oldc]; 34 35 if (newc < 0) { 36 perm[fMax] = oldc; 37 iperm[oldc] = fMax++; 38 } 39 } 40 } 41 if (fMax != fEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of depth %d faces %d does not match permuted nubmer %d", d, fEnd-fStart, fMax-fStart); 42 } 43 *clperm = perm; 44 *invclperm = iperm; 45 PetscFunctionReturn(0); 46 } 47 48 /*@ 49 DMPlexGetOrdering - Calculate a reordering of the mesh 50 51 Collective on dm 52 53 Input Parameter: 54 + dm - The DMPlex object 55 . otype - type of reordering, one of the following: 56 $ MATORDERINGNATURAL - Natural 57 $ MATORDERINGND - Nested Dissection 58 $ MATORDERING1WD - One-way Dissection 59 $ MATORDERINGRCM - Reverse Cuthill-McKee 60 $ MATORDERINGQMD - Quotient Minimum Degree 61 - label - [Optional] Label used to segregate ordering into sets, or NULL 62 63 64 Output Parameter: 65 . perm - The point permutation as an IS, perm[old point number] = new point number 66 67 Note: The label is used to group sets of points together by label value. This makes it easy to reorder a mesh which 68 has different types of cells, and then loop over each set of reordered cells for assembly. 69 70 Level: intermediate 71 72 .seealso: MatGetOrdering() 73 @*/ 74 PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm) 75 { 76 PetscInt numCells = 0; 77 PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i; 78 PetscErrorCode ierr; 79 80 PetscFunctionBegin; 81 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 82 PetscValidPointer(perm, 3); 83 ierr = DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency);CHKERRQ(ierr); 84 ierr = PetscMalloc3(numCells,&cperm,numCells,&mask,numCells*2,&xls);CHKERRQ(ierr); 85 if (numCells) { 86 /* Shift for Fortran numbering */ 87 for (i = 0; i < start[numCells]; ++i) ++adjacency[i]; 88 for (i = 0; i <= numCells; ++i) ++start[i]; 89 ierr = SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls);CHKERRQ(ierr); 90 } 91 ierr = PetscFree(start);CHKERRQ(ierr); 92 ierr = PetscFree(adjacency);CHKERRQ(ierr); 93 /* Shift for Fortran numbering */ 94 for (c = 0; c < numCells; ++c) --cperm[c]; 95 /* Segregate */ 96 if (label) { 97 IS valueIS; 98 const PetscInt *values; 99 PetscInt numValues, numPoints = 0; 100 PetscInt *sperm, *vsize, *voff, v; 101 102 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 103 ierr = ISSort(valueIS);CHKERRQ(ierr); 104 ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr); 105 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 106 ierr = PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff);CHKERRQ(ierr); 107 for (v = 0; v < numValues; ++v) { 108 ierr = DMLabelGetStratumSize(label, values[v], &vsize[v]);CHKERRQ(ierr); 109 if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1]; 110 numPoints += vsize[v]; 111 } 112 if (numPoints != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %D cells < %D total", numPoints, numCells); 113 for (c = 0; c < numCells; ++c) { 114 const PetscInt oldc = cperm[c]; 115 PetscInt val, vloc; 116 117 ierr = DMLabelGetValue(label, oldc, &val);CHKERRQ(ierr); 118 if (val == -1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %D not present in label", oldc); 119 ierr = PetscFindInt(val, numValues, values, &vloc);CHKERRQ(ierr); 120 if (vloc < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D not present label", val); 121 sperm[voff[vloc+1]++] = oldc; 122 } 123 for (v = 0; v < numValues; ++v) { 124 if (voff[v+1] - voff[v] != vsize[v]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of %D values found is %D != %D", values[v], voff[v+1] - voff[v], vsize[v]); 125 } 126 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 127 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 128 ierr = PetscArraycpy(cperm, sperm, numCells);CHKERRQ(ierr); 129 ierr = PetscFree3(sperm, vsize, voff);CHKERRQ(ierr); 130 } 131 /* Construct closure */ 132 ierr = DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm);CHKERRQ(ierr); 133 ierr = PetscFree3(cperm,mask,xls);CHKERRQ(ierr); 134 ierr = PetscFree(clperm);CHKERRQ(ierr); 135 /* Invert permutation */ 136 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 137 ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm);CHKERRQ(ierr); 138 PetscFunctionReturn(0); 139 } 140 141 /*@ 142 DMPlexPermute - Reorder the mesh according to the input permutation 143 144 Collective on dm 145 146 Input Parameter: 147 + dm - The DMPlex object 148 - perm - The point permutation, perm[old point number] = new point number 149 150 Output Parameter: 151 . pdm - The permuted DM 152 153 Level: intermediate 154 155 .seealso: MatPermute() 156 @*/ 157 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) 158 { 159 DM_Plex *plex = (DM_Plex *) dm->data, *plexNew; 160 PetscSection section, sectionNew; 161 PetscInt dim; 162 PetscErrorCode ierr; 163 164 PetscFunctionBegin; 165 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 166 PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 167 PetscValidPointer(pdm, 3); 168 ierr = DMCreate(PetscObjectComm((PetscObject) dm), pdm);CHKERRQ(ierr); 169 ierr = DMSetType(*pdm, DMPLEX);CHKERRQ(ierr); 170 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 171 ierr = DMSetDimension(*pdm, dim);CHKERRQ(ierr); 172 ierr = DMCopyDisc(dm, *pdm);CHKERRQ(ierr); 173 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 174 if (section) { 175 ierr = PetscSectionPermute(section, perm, §ionNew);CHKERRQ(ierr); 176 ierr = DMSetLocalSection(*pdm, sectionNew);CHKERRQ(ierr); 177 ierr = PetscSectionDestroy(§ionNew);CHKERRQ(ierr); 178 } 179 plexNew = (DM_Plex *) (*pdm)->data; 180 /* Ignore ltogmap, ltogmapb */ 181 /* Ignore sf, sectionSF */ 182 /* Ignore globalVertexNumbers, globalCellNumbers */ 183 /* Reorder labels */ 184 { 185 PetscInt numLabels, l; 186 DMLabel label, labelNew; 187 188 ierr = DMGetNumLabels(dm, &numLabels);CHKERRQ(ierr); 189 for (l = 0; l < numLabels; ++l) { 190 ierr = DMGetLabelByNum(dm, l, &label);CHKERRQ(ierr); 191 ierr = DMLabelPermute(label, perm, &labelNew);CHKERRQ(ierr); 192 ierr = DMAddLabel(*pdm, labelNew);CHKERRQ(ierr); 193 ierr = DMLabelDestroy(&labelNew);CHKERRQ(ierr); 194 } 195 ierr = DMGetLabel(*pdm, "depth", &(*pdm)->depthLabel);CHKERRQ(ierr); 196 if (plex->subpointMap) {ierr = DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap);CHKERRQ(ierr);} 197 } 198 /* Reorder topology */ 199 { 200 const PetscInt *pperm; 201 PetscInt maxConeSize, maxSupportSize, n, pStart, pEnd, p; 202 203 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 204 plexNew->maxConeSize = maxConeSize; 205 plexNew->maxSupportSize = maxSupportSize; 206 ierr = PetscSectionDestroy(&plexNew->coneSection);CHKERRQ(ierr); 207 ierr = PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection);CHKERRQ(ierr); 208 ierr = PetscSectionGetStorageSize(plexNew->coneSection, &n);CHKERRQ(ierr); 209 ierr = PetscMalloc1(n, &plexNew->cones);CHKERRQ(ierr); 210 ierr = PetscMalloc1(n, &plexNew->coneOrientations);CHKERRQ(ierr); 211 ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr); 212 ierr = PetscSectionGetChart(plex->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 213 for (p = pStart; p < pEnd; ++p) { 214 PetscInt dof, off, offNew, d; 215 216 ierr = PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof);CHKERRQ(ierr); 217 ierr = PetscSectionGetOffset(plex->coneSection, p, &off);CHKERRQ(ierr); 218 ierr = PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew);CHKERRQ(ierr); 219 for (d = 0; d < dof; ++d) { 220 plexNew->cones[offNew+d] = pperm[plex->cones[off+d]]; 221 plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d]; 222 } 223 } 224 ierr = PetscSectionDestroy(&plexNew->supportSection);CHKERRQ(ierr); 225 ierr = PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection);CHKERRQ(ierr); 226 ierr = PetscSectionGetStorageSize(plexNew->supportSection, &n);CHKERRQ(ierr); 227 ierr = PetscMalloc1(n, &plexNew->supports);CHKERRQ(ierr); 228 ierr = PetscSectionGetChart(plex->supportSection, &pStart, &pEnd);CHKERRQ(ierr); 229 for (p = pStart; p < pEnd; ++p) { 230 PetscInt dof, off, offNew, d; 231 232 ierr = PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof);CHKERRQ(ierr); 233 ierr = PetscSectionGetOffset(plex->supportSection, p, &off);CHKERRQ(ierr); 234 ierr = PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew);CHKERRQ(ierr); 235 for (d = 0; d < dof; ++d) { 236 plexNew->supports[offNew+d] = pperm[plex->supports[off+d]]; 237 } 238 } 239 ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr); 240 } 241 /* Remap coordinates */ 242 { 243 DM cdm, cdmNew; 244 PetscSection csection, csectionNew; 245 Vec coordinates, coordinatesNew; 246 PetscScalar *coords, *coordsNew; 247 const PetscInt *pperm; 248 PetscInt pStart, pEnd, p; 249 const char *name; 250 251 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 252 ierr = DMGetLocalSection(cdm, &csection);CHKERRQ(ierr); 253 ierr = PetscSectionPermute(csection, perm, &csectionNew);CHKERRQ(ierr); 254 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 255 ierr = VecDuplicate(coordinates, &coordinatesNew);CHKERRQ(ierr); 256 ierr = PetscObjectGetName((PetscObject)coordinates,&name);CHKERRQ(ierr); 257 ierr = PetscObjectSetName((PetscObject)coordinatesNew,name);CHKERRQ(ierr); 258 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 259 ierr = VecGetArray(coordinatesNew, &coordsNew);CHKERRQ(ierr); 260 ierr = PetscSectionGetChart(csectionNew, &pStart, &pEnd);CHKERRQ(ierr); 261 ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr); 262 for (p = pStart; p < pEnd; ++p) { 263 PetscInt dof, off, offNew, d; 264 265 ierr = PetscSectionGetDof(csectionNew, p, &dof);CHKERRQ(ierr); 266 ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr); 267 ierr = PetscSectionGetOffset(csectionNew, pperm[p], &offNew);CHKERRQ(ierr); 268 for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d]; 269 } 270 ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr); 271 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 272 ierr = VecRestoreArray(coordinatesNew, &coordsNew);CHKERRQ(ierr); 273 ierr = DMGetCoordinateDM(*pdm, &cdmNew);CHKERRQ(ierr); 274 ierr = DMSetLocalSection(cdmNew, csectionNew);CHKERRQ(ierr); 275 ierr = DMSetCoordinatesLocal(*pdm, coordinatesNew);CHKERRQ(ierr); 276 ierr = PetscSectionDestroy(&csectionNew);CHKERRQ(ierr); 277 ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr); 278 } 279 (*pdm)->setupcalled = PETSC_TRUE; 280 PetscFunctionReturn(0); 281 } 282