1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #ifdef 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);break; 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 = DMGetDefaultSection(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 PetscFunctionList DMPlexGenerateList; 164 165 struct _n_PetscFunctionList { 166 PetscErrorCode (*generate)(DM, PetscBool, DM*); 167 PetscErrorCode (*refine)(DM,double*, DM*); 168 char *name; /* string to identify routine */ 169 PetscInt dim; 170 PetscFunctionList next; /* next pointer */ 171 }; 172 173 PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined) 174 { 175 PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *); 176 PetscReal refinementLimit; 177 PetscInt dim, cStart, cEnd; 178 char genname[1024], *name = NULL; 179 PetscBool flg, localized; 180 PetscErrorCode ierr; 181 PetscErrorCode (*refine)(DM,double*,DM*); 182 PetscFunctionList fl; 183 PetscReal *maxVolumes; 184 PetscInt c; 185 186 PetscFunctionBegin; 187 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 188 ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr); 189 ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr); 190 if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0); 191 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 192 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 193 ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr); 194 if (flg) name = genname; 195 196 fl = DMPlexGenerateList; 197 if (name) { 198 while (fl) { 199 ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr); 200 if (flg) { 201 refine = fl->refine; 202 goto gotit; 203 } 204 fl = fl->next; 205 } 206 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid refiner %g not registered",name);CHKERRQ(ierr); 207 } else { 208 while (fl) { 209 if (dim-1 == fl->dim) { 210 refine = fl->refine; 211 goto gotit; 212 } 213 fl = fl->next; 214 } 215 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid refiner of dimension %D registered",dim);CHKERRQ(ierr); 216 } 217 218 gotit: switch (dim) { 219 case 2: 220 ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 221 if (adaptLabel) { 222 ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr); 223 } else if (refinementFunc) { 224 for (c = cStart; c < cEnd; ++c) { 225 PetscReal vol, centroid[3]; 226 PetscReal maxVol; 227 228 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 229 ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr); 230 maxVolumes[c - cStart] = (double) maxVol; 231 } 232 } else { 233 for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 234 } 235 #if !defined(PETSC_USE_REAL_DOUBLE) 236 { 237 double *mvols; 238 ierr = PetscMalloc1(cEnd - cStart,&mvols);CHKERRQ(ierr); 239 for (c = 0; c < cEnd-cStart; ++c) mvols[c] = (double)maxVolumes[c]; 240 ierr = (*refine)(dm, mvols, dmRefined);CHKERRQ(ierr); 241 ierr = PetscFree(mvols);CHKERRQ(ierr); 242 } 243 #else 244 ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 245 #endif 246 ierr = PetscFree(maxVolumes);CHKERRQ(ierr); 247 break; 248 case 3: 249 ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr); 250 if (adaptLabel) { 251 ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr); 252 } else if (refinementFunc) { 253 for (c = cStart; c < cEnd; ++c) { 254 PetscReal vol, centroid[3]; 255 256 ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr); 257 ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr); 258 } 259 } else { 260 for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit; 261 } 262 ierr = (*refine)(dm, maxVolumes, dmRefined);CHKERRQ(ierr); 263 ierr = PetscFree(maxVolumes);CHKERRQ(ierr); 264 break; 265 default: 266 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim); 267 } 268 ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr); 269 if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);} 270 PetscFunctionReturn(0); 271 } 272 273 PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened) 274 { 275 Vec metricVec; 276 PetscInt cStart, cEnd, vStart, vEnd; 277 DMLabel bdLabel = NULL; 278 char bdLabelName[PETSC_MAX_PATH_LEN]; 279 PetscBool localized, flg; 280 PetscErrorCode ierr; 281 282 PetscFunctionBegin; 283 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 284 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 285 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 286 ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr); 287 ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, PETSC_MAX_PATH_LEN-1, &flg);CHKERRQ(ierr); 288 if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);} 289 ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr); 290 ierr = VecDestroy(&metricVec);CHKERRQ(ierr); 291 if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);} 292 PetscFunctionReturn(0); 293 } 294 295 PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted) 296 { 297 IS flagIS; 298 const PetscInt *flags; 299 PetscInt defFlag, minFlag, maxFlag, numFlags, f; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr); 304 minFlag = defFlag; 305 maxFlag = defFlag; 306 ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr); 307 ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr); 308 ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr); 309 for (f = 0; f < numFlags; ++f) { 310 const PetscInt flag = flags[f]; 311 312 minFlag = PetscMin(minFlag, flag); 313 maxFlag = PetscMax(maxFlag, flag); 314 } 315 ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr); 316 ierr = ISDestroy(&flagIS);CHKERRQ(ierr); 317 { 318 PetscInt minMaxFlag[2], minMaxFlagGlobal[2]; 319 320 minMaxFlag[0] = minFlag; 321 minMaxFlag[1] = -maxFlag; 322 ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 323 minFlag = minMaxFlagGlobal[0]; 324 maxFlag = -minMaxFlagGlobal[1]; 325 } 326 if (minFlag == maxFlag) { 327 switch (minFlag) { 328 case DM_ADAPT_DETERMINE: 329 *dmAdapted = NULL;break; 330 case DM_ADAPT_REFINE: 331 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr); 332 ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break; 333 case DM_ADAPT_COARSEN: 334 ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break; 335 default: 336 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag);break; 337 } 338 } else { 339 ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr); 340 ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr); 341 } 342 PetscFunctionReturn(0); 343 } 344 345 /* 346 DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field. 347 348 Input Parameters: 349 + dm - The DM object 350 . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector 351 - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_". 352 353 Output Parameter: 354 . dmNew - the new DM 355 356 Level: advanced 357 358 .seealso: DMCoarsen(), DMRefine() 359 */ 360 PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew) 361 { 362 #ifdef PETSC_HAVE_PRAGMATIC 363 MPI_Comm comm; 364 const char *bdName = "_boundary_"; 365 #if 0 366 DM odm = dm; 367 #endif 368 DM udm, cdm; 369 DMLabel bdLabelFull; 370 const char *bdLabelName; 371 IS bdIS, globalVertexNum; 372 PetscSection coordSection; 373 Vec coordinates; 374 const PetscScalar *coords, *met; 375 const PetscInt *bdFacesFull, *gV; 376 PetscInt *bdFaces, *bdFaceIds, *l2gv; 377 PetscReal *x, *y, *z, *metric; 378 PetscInt *cells; 379 PetscInt dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v; 380 PetscInt off, maxConeSize, numBdFaces, f, bdSize; 381 PetscBool flg; 382 DMLabel bdLabelNew; 383 double *coordsNew; 384 PetscInt *bdTags; 385 PetscReal *xNew[3] = {NULL, NULL, NULL}; 386 PetscInt *cellsNew; 387 PetscInt d, numCellsNew, numVerticesNew; 388 PetscInt numCornersNew, fStart, fEnd; 389 PetscMPIInt numProcs; 390 PetscErrorCode ierr; 391 392 PetscFunctionBegin; 393 /* Check for FEM adjacency flags */ 394 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 395 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 396 if (bdLabel) { 397 ierr = DMLabelGetName(bdLabel, &bdLabelName);CHKERRQ(ierr); 398 ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr); 399 if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 400 } 401 /* Add overlap for Pragmatic */ 402 #if 0 403 /* Check for overlap by looking for cell in the SF */ 404 if (!overlapped) { 405 ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr); 406 if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);} 407 } 408 #endif 409 /* Get mesh information */ 410 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 411 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 412 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 413 ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 414 ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 415 numCells = cEnd - cStart; 416 numVertices = vEnd - vStart; 417 ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 418 for (c = 0, coff = 0; c < numCells; ++c) { 419 const PetscInt *cone; 420 PetscInt coneSize, cl; 421 422 ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 423 ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 424 for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 425 } 426 ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr); 427 ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr); 428 ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr); 429 for (v = 0, numLocVertices = 0; v < numVertices; ++v) { 430 if (gV[v] >= 0) ++numLocVertices; 431 l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v]; 432 } 433 ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr); 434 ierr = DMDestroy(&udm);CHKERRQ(ierr); 435 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 436 ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr); 437 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 438 ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 439 for (v = vStart; v < vEnd; ++v) { 440 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 441 x[v-vStart] = PetscRealPart(coords[off+0]); 442 if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]); 443 if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]); 444 } 445 ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 446 /* Get boundary mesh */ 447 ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr); 448 ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr); 449 ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr); 450 ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr); 451 ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr); 452 for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 453 PetscInt *closure = NULL; 454 PetscInt closureSize, cl; 455 456 ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 457 for (cl = 0; cl < closureSize*2; cl += 2) { 458 if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 459 } 460 ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 461 } 462 ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 463 for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 464 PetscInt *closure = NULL; 465 PetscInt closureSize, cl; 466 467 ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 468 for (cl = 0; cl < closureSize*2; cl += 2) { 469 if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 470 } 471 ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 472 if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);} 473 else {bdFaceIds[f] = 1;} 474 } 475 ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 476 ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr); 477 /* Get metric */ 478 ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr); 479 ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr); 480 for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]); 481 ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr); 482 #if 0 483 /* Destroy overlap mesh */ 484 ierr = DMDestroy(&dm);CHKERRQ(ierr); 485 #endif 486 /* Create new mesh */ 487 switch (dim) { 488 case 2: 489 pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break; 490 case 3: 491 pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break; 492 default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 493 } 494 pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 495 pragmatic_set_metric(metric); 496 pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0); 497 ierr = PetscFree(l2gv);CHKERRQ(ierr); 498 /* Read out mesh */ 499 pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew); 500 ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr); 501 switch (dim) { 502 case 2: 503 numCornersNew = 3; 504 ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr); 505 pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]); 506 break; 507 case 3: 508 numCornersNew = 4; 509 ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr); 510 pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]); 511 break; 512 default: 513 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 514 } 515 for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = (double) xNew[d][v];} 516 ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr); 517 pragmatic_get_elements(cellsNew); 518 ierr = DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr); 519 /* Read out boundary label */ 520 pragmatic_get_boundaryTags(&bdTags); 521 ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr); 522 ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr); 523 ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr); 524 ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 525 ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr); 526 for (c = cStart; c < cEnd; ++c) { 527 /* Only for simplicial meshes */ 528 coff = (c-cStart)*(dim+1); 529 /* d is the local cell number of the vertex opposite to the face we are marking */ 530 for (d = 0; d < dim+1; ++d) { 531 if (bdTags[coff+d]) { 532 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 */ 533 const PetscInt *cone; 534 535 /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */ 536 ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr); 537 ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr); 538 } 539 } 540 } 541 /* Cleanup */ 542 switch (dim) { 543 case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break; 544 case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break; 545 } 546 ierr = PetscFree(cellsNew);CHKERRQ(ierr); 547 ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 548 ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 549 ierr = PetscFree(coordsNew);CHKERRQ(ierr); 550 pragmatic_finalize(); 551 PetscFunctionReturn(0); 552 #else 553 PetscFunctionBegin; 554 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic."); 555 PetscFunctionReturn(0); 556 #endif 557 } 558