1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #if defined(PETSC_HAVE_PRAGMATIC) 3 #include <pragmatic/cpragmatic.h> 4 #endif 5 6 static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[]) 7 { 8 PetscInt dim, c; 9 PetscErrorCode ierr; 10 11 PetscFunctionBegin; 12 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 13 refRatio = refRatio == PETSC_DEFAULT ? (PetscReal) ((PetscInt) 1 << dim) : refRatio; 14 for (c = cStart; c < cEnd; c++) { 15 PetscReal vol; 16 PetscInt closureSize = 0, cl; 17 PetscInt *closure = NULL; 18 PetscBool anyRefine = PETSC_FALSE; 19 PetscBool anyCoarsen = PETSC_FALSE; 20 PetscBool anyKeep = PETSC_FALSE; 21 22 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr); 23 maxVolumes[c - cStart] = vol; 24 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 25 for (cl = 0; cl < closureSize*2; cl += 2) { 26 const PetscInt point = closure[cl]; 27 PetscInt refFlag; 28 29 ierr = DMLabelGetValue(adaptLabel, point, &refFlag);CHKERRQ(ierr); 30 switch (refFlag) { 31 case DM_ADAPT_REFINE: 32 anyRefine = PETSC_TRUE;break; 33 case DM_ADAPT_COARSEN: 34 anyCoarsen = PETSC_TRUE;break; 35 case DM_ADAPT_KEEP: 36 anyKeep = PETSC_TRUE;break; 37 case DM_ADAPT_DETERMINE: 38 break; 39 default: 40 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %D\n", refFlag); 41 } 42 if (anyRefine) break; 43 } 44 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 45 if (anyRefine) { 46 maxVolumes[c - cStart] = vol / refRatio; 47 } else if (anyKeep) { 48 maxVolumes[c - cStart] = vol; 49 } else if (anyCoarsen) { 50 maxVolumes[c - cStart] = vol * refRatio; 51 } 52 } 53 PetscFunctionReturn(0); 54 } 55 56 static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec) 57 { 58 DM udm, coordDM; 59 PetscSection coordSection; 60 Vec coordinates, mb, mx; 61 Mat A; 62 PetscScalar *metric, *eqns; 63 const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1/refRatio; 64 PetscInt dim, Nv, Neq, c, v; 65 PetscErrorCode ierr; 66 67 PetscFunctionBegin; 68 ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 69 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 70 ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 71 ierr = DMGetLocalSection(coordDM, &coordSection);CHKERRQ(ierr); 72 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 73 Nv = vEnd - vStart; 74 ierr = VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec);CHKERRQ(ierr); 75 ierr = VecGetArray(*metricVec, &metric);CHKERRQ(ierr); 76 Neq = (dim*(dim+1))/2; 77 ierr = PetscMalloc1(PetscSqr(Neq), &eqns);CHKERRQ(ierr); 78 ierr = MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A);CHKERRQ(ierr); 79 ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr); 80 ierr = VecSet(mb, 1.0);CHKERRQ(ierr); 81 for (c = cStart; c < cEnd; ++c) { 82 const PetscScalar *sol; 83 PetscScalar *cellCoords = NULL; 84 PetscReal e[3], vol; 85 const PetscInt *cone; 86 PetscInt coneSize, cl, i, j, d, r; 87 88 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr); 89 /* Only works for simplices */ 90 for (i = 0, r = 0; i < dim+1; ++i) { 91 for (j = 0; j < i; ++j, ++r) { 92 for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i*dim+d] - cellCoords[j*dim+d]); 93 /* FORTRAN ORDERING */ 94 switch (dim) { 95 case 2: 96 eqns[0*Neq+r] = PetscSqr(e[0]); 97 eqns[1*Neq+r] = 2.0*e[0]*e[1]; 98 eqns[2*Neq+r] = PetscSqr(e[1]); 99 break; 100 case 3: 101 eqns[0*Neq+r] = PetscSqr(e[0]); 102 eqns[1*Neq+r] = 2.0*e[0]*e[1]; 103 eqns[2*Neq+r] = 2.0*e[0]*e[2]; 104 eqns[3*Neq+r] = PetscSqr(e[1]); 105 eqns[4*Neq+r] = 2.0*e[1]*e[2]; 106 eqns[5*Neq+r] = PetscSqr(e[2]); 107 break; 108 } 109 } 110 } 111 ierr = MatSetUnfactored(A);CHKERRQ(ierr); 112 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr); 113 ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr); 114 ierr = MatSolve(A, mb, mx);CHKERRQ(ierr); 115 ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr); 116 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr); 117 ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 118 ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 119 for (cl = 0; cl < coneSize; ++cl) { 120 const PetscInt v = cone[cl] - vStart; 121 122 if (dim == 2) { 123 metric[v*4+0] += vol*coarseRatio*sol[0]; 124 metric[v*4+1] += vol*coarseRatio*sol[1]; 125 metric[v*4+2] += vol*coarseRatio*sol[1]; 126 metric[v*4+3] += vol*coarseRatio*sol[2]; 127 } else { 128 metric[v*9+0] += vol*coarseRatio*sol[0]; 129 metric[v*9+1] += vol*coarseRatio*sol[1]; 130 metric[v*9+3] += vol*coarseRatio*sol[1]; 131 metric[v*9+2] += vol*coarseRatio*sol[2]; 132 metric[v*9+6] += vol*coarseRatio*sol[2]; 133 metric[v*9+4] += vol*coarseRatio*sol[3]; 134 metric[v*9+5] += vol*coarseRatio*sol[4]; 135 metric[v*9+7] += vol*coarseRatio*sol[4]; 136 metric[v*9+8] += vol*coarseRatio*sol[5]; 137 } 138 } 139 ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr); 140 } 141 for (v = 0; v < Nv; ++v) { 142 const PetscInt *support; 143 PetscInt supportSize, s; 144 PetscReal vol, totVol = 0.0; 145 146 ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr); 147 ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr); 148 for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;} 149 for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol; 150 } 151 ierr = PetscFree(eqns);CHKERRQ(ierr); 152 ierr = VecRestoreArray(*metricVec, &metric);CHKERRQ(ierr); 153 ierr = VecDestroy(&mx);CHKERRQ(ierr); 154 ierr = VecDestroy(&mb);CHKERRQ(ierr); 155 ierr = MatDestroy(&A);CHKERRQ(ierr); 156 ierr = DMDestroy(&udm);CHKERRQ(ierr); 157 PetscFunctionReturn(0); 158 } 159 160 /* 161 Contains the list of registered DMPlexGenerators routines 162 */ 163 extern PlexGeneratorFunctionList DMPlexGenerateList; 164 165 PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined) 166 { 167 PlexGeneratorFunctionList fl; 168 PetscErrorCode (*refine)(DM,PetscReal*,DM*); 169 PetscErrorCode (*adapt)(DM,DMLabel,DM*); 170 PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *); 171 char genname[PETSC_MAX_PATH_LEN], *name = NULL; 172 PetscReal refinementLimit; 173 PetscReal *maxVolumes; 174 PetscInt dim, cStart, cEnd, c; 175 PetscBool flg, flg2, localized; 176 PetscErrorCode ierr; 177 178 PetscFunctionBegin; 179 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 180 ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr); 181 ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr); 182 if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0); 183 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 184 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 185 ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_adaptor", genname, sizeof(genname), &flg);CHKERRQ(ierr); 186 if (flg) name = genname; 187 else { 188 ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, sizeof(genname), &flg2);CHKERRQ(ierr); 189 if (flg2) name = genname; 190 } 191 192 fl = DMPlexGenerateList; 193 if (name) { 194 while (fl) { 195 ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr); 196 if (flg) { 197 refine = fl->refine; 198 adapt = fl->adaptlabel; 199 goto gotit; 200 } 201 fl = fl->next; 202 } 203 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %s not registered",name); 204 } else { 205 while (fl) { 206 if (fl->dim < 0 || dim-1 == fl->dim) { 207 refine = fl->refine; 208 adapt = fl->adaptlabel; 209 goto gotit; 210 } 211 fl = fl->next; 212 } 213 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %D registered",dim); 214 } 215 216 gotit: 217 switch (dim) { 218 case 1: 219 case 2: 220 case 3: 221 if (adapt) { 222 ierr = (*adapt)(dm, adaptLabel, dmRefined);CHKERRQ(ierr); 223 } else { 224 ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 225 if (adaptLabel) { 226 ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr); 227 } else if (refinementFunc) { 228 for (c = cStart; c < cEnd; ++c) { 229 PetscReal vol, centroid[3]; 230 231 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 232 ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr); 233 } 234 } else { 235 for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 236 } 237 ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 238 ierr = PetscFree(maxVolumes);CHKERRQ(ierr); 239 } 240 break; 241 default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %D is not supported.", dim); 242 } 243 ((DM_Plex *) (*dmRefined)->data)->useHashLocation = ((DM_Plex *) dm->data)->useHashLocation; 244 ierr = DMCopyDisc(dm, *dmRefined);CHKERRQ(ierr); 245 if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);} 246 PetscFunctionReturn(0); 247 } 248 249 PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened) 250 { 251 Vec metricVec; 252 PetscInt cStart, cEnd, vStart, vEnd; 253 DMLabel bdLabel = NULL; 254 char bdLabelName[PETSC_MAX_PATH_LEN]; 255 PetscBool localized, flg; 256 PetscErrorCode ierr; 257 258 PetscFunctionBegin; 259 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 260 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 261 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 262 ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr); 263 ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, sizeof(bdLabelName), &flg);CHKERRQ(ierr); 264 if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);} 265 ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr); 266 ierr = VecDestroy(&metricVec);CHKERRQ(ierr); 267 ((DM_Plex *) (*dmCoarsened)->data)->useHashLocation = ((DM_Plex *) dm->data)->useHashLocation; 268 ierr = DMCopyDisc(dm, *dmCoarsened);CHKERRQ(ierr); 269 if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);} 270 PetscFunctionReturn(0); 271 } 272 273 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted) 274 { 275 IS flagIS; 276 const PetscInt *flags; 277 PetscInt defFlag, minFlag, maxFlag, numFlags, f; 278 PetscErrorCode ierr; 279 280 PetscFunctionBegin; 281 ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr); 282 minFlag = defFlag; 283 maxFlag = defFlag; 284 ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr); 285 ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr); 286 ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr); 287 for (f = 0; f < numFlags; ++f) { 288 const PetscInt flag = flags[f]; 289 290 minFlag = PetscMin(minFlag, flag); 291 maxFlag = PetscMax(maxFlag, flag); 292 } 293 ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr); 294 ierr = ISDestroy(&flagIS);CHKERRQ(ierr); 295 { 296 PetscInt minMaxFlag[2], minMaxFlagGlobal[2]; 297 298 minMaxFlag[0] = minFlag; 299 minMaxFlag[1] = -maxFlag; 300 ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr); 301 minFlag = minMaxFlagGlobal[0]; 302 maxFlag = -minMaxFlagGlobal[1]; 303 } 304 if (minFlag == maxFlag) { 305 switch (minFlag) { 306 case DM_ADAPT_DETERMINE: 307 *dmAdapted = NULL;break; 308 case DM_ADAPT_REFINE: 309 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr); 310 ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break; 311 case DM_ADAPT_COARSEN: 312 ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break; 313 default: 314 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag); 315 } 316 } else { 317 ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr); 318 ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr); 319 } 320 PetscFunctionReturn(0); 321 } 322 323 /* 324 DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field. 325 326 Input Parameters: 327 + dm - The DM object 328 . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector 329 - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_". 330 331 Output Parameter: 332 . dmNew - the new DM 333 334 Level: advanced 335 336 .seealso: DMCoarsen(), DMRefine() 337 */ 338 PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew) 339 { 340 #if defined(PETSC_HAVE_PRAGMATIC) 341 MPI_Comm comm; 342 const char *bdName = "_boundary_"; 343 #if 0 344 DM odm = dm; 345 #endif 346 DM udm, cdm; 347 DMLabel bdLabelFull; 348 const char *bdLabelName; 349 IS bdIS, globalVertexNum; 350 PetscSection coordSection; 351 Vec coordinates; 352 const PetscScalar *coords, *met; 353 const PetscInt *bdFacesFull, *gV; 354 PetscInt *bdFaces, *bdFaceIds, *l2gv; 355 PetscReal *x, *y, *z, *metric; 356 PetscInt *cells; 357 PetscInt dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v; 358 PetscInt off, maxConeSize, numBdFaces, f, bdSize; 359 PetscBool flg; 360 DMLabel bdLabelNew; 361 PetscReal *coordsNew; 362 PetscInt *bdTags; 363 PetscReal *xNew[3] = {NULL, NULL, NULL}; 364 PetscInt *cellsNew; 365 PetscInt d, numCellsNew, numVerticesNew; 366 PetscInt numCornersNew, fStart, fEnd; 367 PetscMPIInt numProcs; 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 /* Check for FEM adjacency flags */ 372 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 373 ierr = MPI_Comm_size(comm, &numProcs);CHKERRMPI(ierr); 374 if (bdLabel) { 375 ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr); 376 ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr); 377 if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 378 } 379 /* Add overlap for Pragmatic */ 380 #if 0 381 /* Check for overlap by looking for cell in the SF */ 382 if (!overlapped) { 383 ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr); 384 if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);} 385 } 386 #endif 387 /* Get mesh information */ 388 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 389 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 390 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 391 ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 392 ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 393 numCells = cEnd - cStart; 394 if (numCells == 0) { 395 PetscMPIInt rank; 396 397 ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 398 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot perform mesh adaptation because process %d does not own any cells.", rank); 399 } 400 numVertices = vEnd - vStart; 401 ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 402 for (c = 0, coff = 0; c < numCells; ++c) { 403 const PetscInt *cone; 404 PetscInt coneSize, cl; 405 406 ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 407 ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 408 for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 409 } 410 ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr); 411 ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr); 412 ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr); 413 for (v = 0, numLocVertices = 0; v < numVertices; ++v) { 414 if (gV[v] >= 0) ++numLocVertices; 415 l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v]; 416 } 417 ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr); 418 ierr = DMDestroy(&udm);CHKERRQ(ierr); 419 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 420 ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 421 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 422 ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 423 for (v = vStart; v < vEnd; ++v) { 424 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 425 x[v-vStart] = PetscRealPart(coords[off+0]); 426 if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]); 427 if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]); 428 } 429 ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 430 /* Get boundary mesh */ 431 ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull);CHKERRQ(ierr); 432 ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr); 433 ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr); 434 ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr); 435 ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr); 436 for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 437 PetscInt *closure = NULL; 438 PetscInt closureSize, cl; 439 440 ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 441 for (cl = 0; cl < closureSize*2; cl += 2) { 442 if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 443 } 444 ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 445 } 446 ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 447 for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 448 PetscInt *closure = NULL; 449 PetscInt closureSize, cl; 450 451 ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 452 for (cl = 0; cl < closureSize*2; cl += 2) { 453 if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 454 } 455 ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 456 if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);} 457 else {bdFaceIds[f] = 1;} 458 } 459 ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 460 ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr); 461 /* Get metric */ 462 ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr); 463 ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr); 464 for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]); 465 ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr); 466 #if 0 467 /* Destroy overlap mesh */ 468 ierr = DMDestroy(&dm);CHKERRQ(ierr); 469 #endif 470 /* Create new mesh */ 471 switch (dim) { 472 case 2: 473 pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break; 474 case 3: 475 pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break; 476 default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim); 477 } 478 pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 479 pragmatic_set_metric(metric); 480 pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0); 481 ierr = PetscFree(l2gv);CHKERRQ(ierr); 482 /* Read out mesh */ 483 pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew); 484 ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr); 485 switch (dim) { 486 case 2: 487 numCornersNew = 3; 488 ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr); 489 pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]); 490 break; 491 case 3: 492 numCornersNew = 4; 493 ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr); 494 pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]); 495 break; 496 default: 497 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim); 498 } 499 for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];} 500 ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr); 501 pragmatic_get_elements(cellsNew); 502 ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr); 503 /* Read out boundary label */ 504 pragmatic_get_boundaryTags(&bdTags); 505 ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr); 506 ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr); 507 ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr); 508 ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 509 ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr); 510 for (c = cStart; c < cEnd; ++c) { 511 /* Only for simplicial meshes */ 512 coff = (c-cStart)*(dim+1); 513 /* d is the local cell number of the vertex opposite to the face we are marking */ 514 for (d = 0; d < dim+1; ++d) { 515 if (bdTags[coff+d]) { 516 const PetscInt perm[4][4] = {{-1, -1, -1, -1}, {-1, -1, -1, -1}, {1, 2, 0, -1}, {3, 2, 1, 0}}; /* perm[d] = face opposite */ 517 const PetscInt *cone; 518 519 /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */ 520 ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr); 521 ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr); 522 } 523 } 524 } 525 /* Cleanup */ 526 switch (dim) { 527 case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break; 528 case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break; 529 } 530 ierr = PetscFree(cellsNew);CHKERRQ(ierr); 531 ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 532 ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 533 ierr = PetscFree(coordsNew);CHKERRQ(ierr); 534 pragmatic_finalize(); 535 PetscFunctionReturn(0); 536 #else 537 PetscFunctionBegin; 538 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic."); 539 #endif 540 } 541