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 .keywords: mesh 73 .seealso: MatGetOrdering() 74 @*/ 75 PetscErrorCode DMPlexGetOrdering(DM dm, MatOrderingType otype, DMLabel label, IS *perm) 76 { 77 PetscInt numCells = 0; 78 PetscInt *start = NULL, *adjacency = NULL, *cperm, *clperm = NULL, *invclperm = NULL, *mask, *xls, pStart, pEnd, c, i; 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 83 PetscValidPointer(perm, 3); 84 ierr = DMPlexCreateNeighborCSR(dm, 0, &numCells, &start, &adjacency);CHKERRQ(ierr); 85 ierr = PetscMalloc3(numCells,&cperm,numCells,&mask,numCells*2,&xls);CHKERRQ(ierr); 86 if (numCells) { 87 /* Shift for Fortran numbering */ 88 for (i = 0; i < start[numCells]; ++i) ++adjacency[i]; 89 for (i = 0; i <= numCells; ++i) ++start[i]; 90 ierr = SPARSEPACKgenrcm(&numCells, start, adjacency, cperm, mask, xls);CHKERRQ(ierr); 91 } 92 ierr = PetscFree(start);CHKERRQ(ierr); 93 ierr = PetscFree(adjacency);CHKERRQ(ierr); 94 /* Shift for Fortran numbering */ 95 for (c = 0; c < numCells; ++c) --cperm[c]; 96 /* Segregate */ 97 if (label) { 98 IS valueIS; 99 const PetscInt *values; 100 PetscInt numValues, numPoints = 0; 101 PetscInt *sperm, *vsize, *voff, v; 102 103 ierr = DMLabelGetValueIS(label, &valueIS);CHKERRQ(ierr); 104 ierr = ISSort(valueIS);CHKERRQ(ierr); 105 ierr = ISGetLocalSize(valueIS, &numValues);CHKERRQ(ierr); 106 ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 107 ierr = PetscCalloc3(numCells,&sperm,numValues,&vsize,numValues+1,&voff);CHKERRQ(ierr); 108 for (v = 0; v < numValues; ++v) { 109 ierr = DMLabelGetStratumSize(label, values[v], &vsize[v]);CHKERRQ(ierr); 110 if (v < numValues-1) voff[v+2] += vsize[v] + voff[v+1]; 111 numPoints += vsize[v]; 112 } 113 if (numPoints != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label only covers %D cells < %D total", numPoints, numCells); 114 for (c = 0; c < numCells; ++c) { 115 const PetscInt oldc = cperm[c]; 116 PetscInt val, vloc; 117 118 ierr = DMLabelGetValue(label, oldc, &val);CHKERRQ(ierr); 119 if (val == -1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %D not present in label", oldc); 120 ierr = PetscFindInt(val, numValues, values, &vloc);CHKERRQ(ierr); 121 if (vloc < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D not present label", val); 122 sperm[voff[vloc+1]++] = oldc; 123 } 124 for (v = 0; v < numValues; ++v) { 125 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]); 126 } 127 ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 128 ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 129 ierr = PetscMemcpy(cperm, sperm, numCells * sizeof(PetscInt));CHKERRQ(ierr); 130 ierr = PetscFree3(sperm, vsize, voff);CHKERRQ(ierr); 131 } 132 /* Construct closure */ 133 ierr = DMPlexCreateOrderingClosure_Static(dm, numCells, cperm, &clperm, &invclperm);CHKERRQ(ierr); 134 ierr = PetscFree3(cperm,mask,xls);CHKERRQ(ierr); 135 ierr = PetscFree(clperm);CHKERRQ(ierr); 136 /* Invert permutation */ 137 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 138 ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd-pStart, invclperm, PETSC_OWN_POINTER, perm);CHKERRQ(ierr); 139 PetscFunctionReturn(0); 140 } 141 142 /*@ 143 DMPlexPermute - Reorder the mesh according to the input permutation 144 145 Collective on DM 146 147 Input Parameter: 148 + dm - The DMPlex object 149 - perm - The point permutation, perm[old point number] = new point number 150 151 Output Parameter: 152 . pdm - The permuted DM 153 154 Level: intermediate 155 156 .keywords: mesh 157 .seealso: MatPermute() 158 @*/ 159 PetscErrorCode DMPlexPermute(DM dm, IS perm, DM *pdm) 160 { 161 DM_Plex *plex = (DM_Plex *) dm->data, *plexNew; 162 PetscSection section, sectionNew; 163 PetscInt dim; 164 PetscErrorCode ierr; 165 166 PetscFunctionBegin; 167 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 168 PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 169 PetscValidPointer(pdm, 3); 170 ierr = DMCreate(PetscObjectComm((PetscObject) dm), pdm);CHKERRQ(ierr); 171 ierr = DMSetType(*pdm, DMPLEX);CHKERRQ(ierr); 172 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 173 ierr = DMSetDimension(*pdm, dim);CHKERRQ(ierr); 174 ierr = DMCopyDisc(dm, *pdm);CHKERRQ(ierr); 175 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 176 if (section) { 177 ierr = PetscSectionPermute(section, perm, §ionNew);CHKERRQ(ierr); 178 ierr = DMSetSection(*pdm, sectionNew);CHKERRQ(ierr); 179 ierr = PetscSectionDestroy(§ionNew);CHKERRQ(ierr); 180 } 181 plexNew = (DM_Plex *) (*pdm)->data; 182 /* Ignore ltogmap, ltogmapb */ 183 /* Ignore sf, defaultSF */ 184 /* Ignore globalVertexNumbers, globalCellNumbers */ 185 /* Remap coordinates */ 186 { 187 DM cdm, cdmNew; 188 PetscSection csection, csectionNew; 189 Vec coordinates, coordinatesNew; 190 PetscScalar *coords, *coordsNew; 191 const PetscInt *pperm; 192 PetscInt pStart, pEnd, p; 193 const char *name; 194 195 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 196 ierr = DMGetSection(cdm, &csection);CHKERRQ(ierr); 197 ierr = PetscSectionPermute(csection, perm, &csectionNew);CHKERRQ(ierr); 198 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 199 ierr = VecDuplicate(coordinates, &coordinatesNew);CHKERRQ(ierr); 200 ierr = PetscObjectGetName((PetscObject)coordinates,&name);CHKERRQ(ierr); 201 ierr = PetscObjectSetName((PetscObject)coordinatesNew,name);CHKERRQ(ierr); 202 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 203 ierr = VecGetArray(coordinatesNew, &coordsNew);CHKERRQ(ierr); 204 ierr = PetscSectionGetChart(csectionNew, &pStart, &pEnd);CHKERRQ(ierr); 205 ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr); 206 for (p = pStart; p < pEnd; ++p) { 207 PetscInt dof, off, offNew, d; 208 209 ierr = PetscSectionGetDof(csectionNew, p, &dof);CHKERRQ(ierr); 210 ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr); 211 ierr = PetscSectionGetOffset(csectionNew, pperm[p], &offNew);CHKERRQ(ierr); 212 for (d = 0; d < dof; ++d) coordsNew[offNew+d] = coords[off+d]; 213 } 214 ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr); 215 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 216 ierr = VecRestoreArray(coordinatesNew, &coordsNew);CHKERRQ(ierr); 217 ierr = DMGetCoordinateDM(*pdm, &cdmNew);CHKERRQ(ierr); 218 ierr = DMSetSection(cdmNew, csectionNew);CHKERRQ(ierr); 219 ierr = DMSetCoordinatesLocal(*pdm, coordinatesNew);CHKERRQ(ierr); 220 ierr = PetscSectionDestroy(&csectionNew);CHKERRQ(ierr); 221 ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr); 222 } 223 /* Reorder labels */ 224 { 225 PetscInt numLabels, l; 226 DMLabel label, labelNew; 227 228 ierr = DMGetNumLabels(dm, &numLabels);CHKERRQ(ierr); 229 for (l = numLabels-1; l >= 0; --l) { 230 ierr = DMGetLabelByNum(dm, l, &label);CHKERRQ(ierr); 231 ierr = DMLabelPermute(label, perm, &labelNew);CHKERRQ(ierr); 232 ierr = DMAddLabel(*pdm, labelNew);CHKERRQ(ierr); 233 } 234 if (plex->subpointMap) {ierr = DMLabelPermute(plex->subpointMap, perm, &plexNew->subpointMap);CHKERRQ(ierr);} 235 } 236 /* Reorder topology */ 237 { 238 const PetscInt *pperm; 239 PetscInt maxConeSize, maxSupportSize, n, pStart, pEnd, p; 240 241 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 242 plexNew->maxConeSize = maxConeSize; 243 plexNew->maxSupportSize = maxSupportSize; 244 ierr = PetscSectionDestroy(&plexNew->coneSection);CHKERRQ(ierr); 245 ierr = PetscSectionPermute(plex->coneSection, perm, &plexNew->coneSection);CHKERRQ(ierr); 246 ierr = PetscSectionGetStorageSize(plexNew->coneSection, &n);CHKERRQ(ierr); 247 ierr = PetscMalloc1(n, &plexNew->cones);CHKERRQ(ierr); 248 ierr = PetscMalloc1(n, &plexNew->coneOrientations);CHKERRQ(ierr); 249 ierr = ISGetIndices(perm, &pperm);CHKERRQ(ierr); 250 ierr = PetscSectionGetChart(plex->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 251 for (p = pStart; p < pEnd; ++p) { 252 PetscInt dof, off, offNew, d; 253 254 ierr = PetscSectionGetDof(plexNew->coneSection, pperm[p], &dof);CHKERRQ(ierr); 255 ierr = PetscSectionGetOffset(plex->coneSection, p, &off);CHKERRQ(ierr); 256 ierr = PetscSectionGetOffset(plexNew->coneSection, pperm[p], &offNew);CHKERRQ(ierr); 257 for (d = 0; d < dof; ++d) { 258 plexNew->cones[offNew+d] = pperm[plex->cones[off+d]]; 259 plexNew->coneOrientations[offNew+d] = plex->coneOrientations[off+d]; 260 } 261 } 262 ierr = PetscSectionDestroy(&plexNew->supportSection);CHKERRQ(ierr); 263 ierr = PetscSectionPermute(plex->supportSection, perm, &plexNew->supportSection);CHKERRQ(ierr); 264 ierr = PetscSectionGetStorageSize(plexNew->supportSection, &n);CHKERRQ(ierr); 265 ierr = PetscMalloc1(n, &plexNew->supports);CHKERRQ(ierr); 266 ierr = PetscSectionGetChart(plex->supportSection, &pStart, &pEnd);CHKERRQ(ierr); 267 for (p = pStart; p < pEnd; ++p) { 268 PetscInt dof, off, offNew, d; 269 270 ierr = PetscSectionGetDof(plexNew->supportSection, pperm[p], &dof);CHKERRQ(ierr); 271 ierr = PetscSectionGetOffset(plex->supportSection, p, &off);CHKERRQ(ierr); 272 ierr = PetscSectionGetOffset(plexNew->supportSection, pperm[p], &offNew);CHKERRQ(ierr); 273 for (d = 0; d < dof; ++d) { 274 plexNew->supports[offNew+d] = pperm[plex->supports[off+d]]; 275 } 276 } 277 ierr = ISRestoreIndices(perm, &pperm);CHKERRQ(ierr); 278 } 279 ierr = DMCopyDisc(dm, *pdm);CHKERRQ(ierr); 280 (*pdm)->setupcalled = PETSC_TRUE; 281 PetscFunctionReturn(0); 282 } 283