Actual source code: plextransform.c

  1: #include <petsc/private/dmplextransformimpl.h>

  3: #include <petsc/private/petscfeimpl.h>

  5: PetscClassId DMPLEXTRANSFORM_CLASSID;

  7: PetscFunctionList DMPlexTransformList              = NULL;
  8: PetscBool         DMPlexTransformRegisterAllCalled = PETSC_FALSE;

 10: PetscLogEvent DMPLEXTRANSFORM_SetUp, DMPLEXTRANSFORM_Apply, DMPLEXTRANSFORM_SetConeSizes, DMPLEXTRANSFORM_SetCones, DMPLEXTRANSFORM_CreateSF, DMPLEXTRANSFORM_CreateLabels, DMPLEXTRANSFORM_SetCoordinates;

 12: /* Construct cell type order since we must loop over cell types in the same dimensional order they are stored in the plex if dm != NULL
 13:         OR in standard plex ordering if dm == NULL */
 14: static PetscErrorCode DMPlexCreateCellTypeOrder_Internal(DM dm, PetscInt dim, PetscInt *ctOrder[], PetscInt *ctOrderInv[])
 15: {
 16:   PetscInt *ctO, *ctOInv;
 17:   PetscInt  d, c, off = 0;
 18:   PetscInt  dimOrder[5] = {3, 2, 1, 0, -1};

 20:   PetscFunctionBegin;
 21:   PetscCall(PetscCalloc2(DM_NUM_POLYTOPES + 1, &ctO, DM_NUM_POLYTOPES + 1, &ctOInv));
 22:   if (dm) { // Order the dimensions by their starting location
 23:     PetscInt hStart[4] = {-1, -1, -1, -1};
 24:     for (d = 0; d <= dim; ++d) PetscCall(DMPlexGetDepthStratum(dm, dim - d, &hStart[d], NULL));
 25:     PetscCall(PetscSortIntWithArray(dim + 1, hStart, &dimOrder[3 - dim]));
 26:   } else if (dim > 1) { // Standard plex ordering. dimOrder is in correct order if dim > 1
 27:     off             = 4 - dim;
 28:     dimOrder[off++] = 0;
 29:     for (d = dim - 1; d > 0; --d) dimOrder[off++] = d;
 30:   }

 32:   off = 0;
 33:   for (d = 0; d < 5; ++d) {
 34:     for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
 35:       if (c == DM_POLYTOPE_UNKNOWN_CELL || c == DM_POLYTOPE_UNKNOWN_FACE) continue;
 36:       if (DMPolytopeTypeGetDim((DMPolytopeType)c) == dimOrder[d]) ctO[off++] = c;
 37:     }
 38:   }
 39:   for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
 40:     if (c == DM_POLYTOPE_UNKNOWN_CELL || c == DM_POLYTOPE_UNKNOWN_FACE) ctO[off++] = c;
 41:   }
 42:   ctO[off++] = DM_NUM_POLYTOPES;
 43:   PetscCheck(off == DM_NUM_POLYTOPES + 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid offset %" PetscInt_FMT " for cell type order", off);

 45:   for (c = 0; c <= DM_NUM_POLYTOPES; ++c) ctOInv[ctO[c]] = c;

 47:   *ctOrder    = ctO;
 48:   *ctOrderInv = ctOInv;
 49:   PetscFunctionReturn(PETSC_SUCCESS);
 50: }

 52: /*@C
 53:   DMPlexTransformRegister - Adds a new transform component implementation

 55:   Not Collective

 57:   Input Parameters:
 58: + name        - The name of a new user-defined creation routine
 59: - create_func - The creation routine

 61:   Example Usage:
 62: .vb
 63:   DMPlexTransformRegister("my_transform", MyTransformCreate);
 64: .ve

 66:   Then, your transform type can be chosen with the procedural interface via
 67: .vb
 68:   DMPlexTransformCreate(MPI_Comm, DMPlexTransform *);
 69:   DMPlexTransformSetType(DMPlexTransform, "my_transform");
 70: .ve
 71:   or at runtime via the option
 72: .vb
 73:   -dm_plex_transform_type my_transform
 74: .ve

 76:   Level: advanced

 78:   Note:
 79:   `DMPlexTransformRegister()` may be called multiple times to add several user-defined transforms

 81: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformRegisterAll()`, `DMPlexTransformRegisterDestroy()`
 82: @*/
 83: PetscErrorCode DMPlexTransformRegister(const char name[], PetscErrorCode (*create_func)(DMPlexTransform))
 84: {
 85:   PetscFunctionBegin;
 86:   PetscCall(DMInitializePackage());
 87:   PetscCall(PetscFunctionListAdd(&DMPlexTransformList, name, create_func));
 88:   PetscFunctionReturn(PETSC_SUCCESS);
 89: }

 91: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Filter(DMPlexTransform);
 92: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Regular(DMPlexTransform);
 93: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_ToBox(DMPlexTransform);
 94: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_ToSimplex(DMPlexTransform);
 95: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Alfeld(DMPlexTransform);
 96: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_SBR(DMPlexTransform);
 97: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_BL(DMPlexTransform);
 98: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_1D(DMPlexTransform);
 99: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Extrude(DMPlexTransform);
100: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Cohesive(DMPlexTransform);

102: /*@C
103:   DMPlexTransformRegisterAll - Registers all of the transform components in the `DM` package.

105:   Not Collective

107:   Level: advanced

109: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRegisterAll()`, `DMPlexTransformRegisterDestroy()`
110: @*/
111: PetscErrorCode DMPlexTransformRegisterAll(void)
112: {
113:   PetscFunctionBegin;
114:   if (DMPlexTransformRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
115:   DMPlexTransformRegisterAllCalled = PETSC_TRUE;

117:   PetscCall(DMPlexTransformRegister(DMPLEXTRANSFORMFILTER, DMPlexTransformCreate_Filter));
118:   PetscCall(DMPlexTransformRegister(DMPLEXREFINEREGULAR, DMPlexTransformCreate_Regular));
119:   PetscCall(DMPlexTransformRegister(DMPLEXREFINETOBOX, DMPlexTransformCreate_ToBox));
120:   PetscCall(DMPlexTransformRegister(DMPLEXREFINETOSIMPLEX, DMPlexTransformCreate_ToSimplex));
121:   PetscCall(DMPlexTransformRegister(DMPLEXREFINEALFELD, DMPlexTransformCreate_Alfeld));
122:   PetscCall(DMPlexTransformRegister(DMPLEXREFINEBOUNDARYLAYER, DMPlexTransformCreate_BL));
123:   PetscCall(DMPlexTransformRegister(DMPLEXREFINESBR, DMPlexTransformCreate_SBR));
124:   PetscCall(DMPlexTransformRegister(DMPLEXREFINE1D, DMPlexTransformCreate_1D));
125:   PetscCall(DMPlexTransformRegister(DMPLEXEXTRUDETYPE, DMPlexTransformCreate_Extrude));
126:   PetscCall(DMPlexTransformRegister(DMPLEXCOHESIVEEXTRUDE, DMPlexTransformCreate_Cohesive));
127:   PetscFunctionReturn(PETSC_SUCCESS);
128: }

130: /*@C
131:   DMPlexTransformRegisterDestroy - This function destroys the registered `DMPlexTransformType`. It is called from `PetscFinalize()`.

133:   Not collective

135:   Level: developer

137: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRegisterAll()`, `DMPlexTransformType`, `PetscInitialize()`
138: @*/
139: PetscErrorCode DMPlexTransformRegisterDestroy(void)
140: {
141:   PetscFunctionBegin;
142:   PetscCall(PetscFunctionListDestroy(&DMPlexTransformList));
143:   DMPlexTransformRegisterAllCalled = PETSC_FALSE;
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }

147: /*@
148:   DMPlexTransformCreate - Creates an empty transform object. The type can then be set with `DMPlexTransformSetType()`.

150:   Collective

152:   Input Parameter:
153: . comm - The communicator for the transform object

155:   Output Parameter:
156: . tr - The transform object

158:   Level: beginner

160: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformType`, `DMPlexTransformSetType()`, `DMPLEXREFINEREGULAR`, `DMPLEXTRANSFORMFILTER`
161: @*/
162: PetscErrorCode DMPlexTransformCreate(MPI_Comm comm, DMPlexTransform *tr)
163: {
164:   DMPlexTransform t;

166:   PetscFunctionBegin;
167:   PetscAssertPointer(tr, 2);
168:   *tr = NULL;
169:   PetscCall(DMInitializePackage());

171:   PetscCall(PetscHeaderCreate(t, DMPLEXTRANSFORM_CLASSID, "DMPlexTransform", "Mesh Transform", "DMPlexTransform", comm, DMPlexTransformDestroy, DMPlexTransformView));
172:   t->setupcalled = PETSC_FALSE;
173:   t->redFactor   = 2.0;
174:   PetscCall(PetscCalloc2(DM_NUM_POLYTOPES, &t->coordFE, DM_NUM_POLYTOPES, &t->refGeom));
175:   *tr = t;
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: /*@
180:   DMPlexTransformSetType - Sets the particular implementation for a transform.

182:   Collective

184:   Input Parameters:
185: + tr     - The transform
186: - method - The name of the transform type

188:   Options Database Key:
189: . -dm_plex_transform_type type - Sets the transform type; see `DMPlexTransformType`

191:   Level: intermediate

193: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformType`, `DMPlexTransformGetType()`, `DMPlexTransformCreate()`
194: @*/
195: PetscErrorCode DMPlexTransformSetType(DMPlexTransform tr, DMPlexTransformType method)
196: {
197:   PetscErrorCode (*r)(DMPlexTransform);
198:   PetscBool match;

200:   PetscFunctionBegin;
202:   PetscCall(PetscObjectTypeCompare((PetscObject)tr, method, &match));
203:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

205:   PetscCall(DMPlexTransformRegisterAll());
206:   PetscCall(PetscFunctionListFind(DMPlexTransformList, method, &r));
207:   PetscCheck(r, PetscObjectComm((PetscObject)tr), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DMPlexTransform type: %s", method);

209:   PetscTryTypeMethod(tr, destroy);
210:   PetscCall(PetscMemzero(tr->ops, sizeof(*tr->ops)));
211:   PetscCall(PetscObjectChangeTypeName((PetscObject)tr, method));
212:   PetscCall((*r)(tr));
213:   PetscFunctionReturn(PETSC_SUCCESS);
214: }

216: /*@
217:   DMPlexTransformGetType - Gets the type name (as a string) from the transform.

219:   Not Collective

221:   Input Parameter:
222: . tr - The `DMPlexTransform`

224:   Output Parameter:
225: . type - The `DMPlexTransformType` name

227:   Level: intermediate

229: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformType`, `DMPlexTransformSetType()`, `DMPlexTransformCreate()`
230: @*/
231: PetscErrorCode DMPlexTransformGetType(DMPlexTransform tr, DMPlexTransformType *type)
232: {
233:   PetscFunctionBegin;
235:   PetscAssertPointer(type, 2);
236:   PetscCall(DMPlexTransformRegisterAll());
237:   *type = ((PetscObject)tr)->type_name;
238:   PetscFunctionReturn(PETSC_SUCCESS);
239: }

241: static PetscErrorCode DMPlexTransformView_Ascii(DMPlexTransform tr, PetscViewer v)
242: {
243:   PetscViewerFormat format;

245:   PetscFunctionBegin;
246:   PetscCall(PetscViewerGetFormat(v, &format));
247:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
248:     const PetscInt *trTypes = NULL;
249:     IS              trIS;
250:     PetscInt        cols = 8;
251:     PetscInt        Nrt  = 8, f, g;

253:     if (tr->trType) PetscCall(DMLabelView(tr->trType, v));
254:     PetscCall(PetscViewerASCIIPrintf(v, "Source Starts\n"));
255:     for (g = 0; g <= cols; ++g) PetscCall(PetscViewerASCIIPrintf(v, " %14s", DMPolytopeTypes[g]));
256:     PetscCall(PetscViewerASCIIPrintf(v, "\n"));
257:     for (f = 0; f <= cols; ++f) PetscCall(PetscViewerASCIIPrintf(v, " %14" PetscInt_FMT, tr->ctStart[f]));
258:     PetscCall(PetscViewerASCIIPrintf(v, "\n"));
259:     PetscCall(PetscViewerASCIIPrintf(v, "Target Starts\n"));
260:     for (g = 0; g <= cols; ++g) PetscCall(PetscViewerASCIIPrintf(v, " %14s", DMPolytopeTypes[g]));
261:     PetscCall(PetscViewerASCIIPrintf(v, "\n"));
262:     for (f = 0; f <= cols; ++f) PetscCall(PetscViewerASCIIPrintf(v, " %14" PetscInt_FMT, tr->ctStartNew[f]));
263:     PetscCall(PetscViewerASCIIPrintf(v, "\n"));

265:     if (tr->trType) {
266:       PetscCall(DMLabelGetNumValues(tr->trType, &Nrt));
267:       PetscCall(DMLabelGetValueIS(tr->trType, &trIS));
268:       PetscCall(ISGetIndices(trIS, &trTypes));
269:     }
270:     PetscCall(PetscViewerASCIIPrintf(v, "Offsets\n"));
271:     PetscCall(PetscViewerASCIIPrintf(v, "     "));
272:     for (g = 0; g < cols; ++g) PetscCall(PetscViewerASCIIPrintf(v, " %14s", DMPolytopeTypes[g]));
273:     PetscCall(PetscViewerASCIIPrintf(v, "\n"));
274:     for (f = 0; f < Nrt; ++f) {
275:       PetscCall(PetscViewerASCIIPrintf(v, "%2" PetscInt_FMT "  |", trTypes ? trTypes[f] : f));
276:       for (g = 0; g < cols; ++g) PetscCall(PetscViewerASCIIPrintf(v, " %14" PetscInt_FMT, tr->offset[f * DM_NUM_POLYTOPES + g]));
277:       PetscCall(PetscViewerASCIIPrintf(v, " |\n"));
278:     }
279:     if (tr->trType) {
280:       PetscCall(ISRestoreIndices(trIS, &trTypes));
281:       PetscCall(ISDestroy(&trIS));
282:     }
283:   }
284:   PetscFunctionReturn(PETSC_SUCCESS);
285: }

287: /*@
288:   DMPlexTransformView - Views a `DMPlexTransform`

290:   Collective

292:   Input Parameters:
293: + tr - the `DMPlexTransform` object to view
294: - v  - the viewer

296:   Level: beginner

298: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformType`, `PetscViewer`, `DMPlexTransformDestroy()`, `DMPlexTransformCreate()`
299: @*/
300: PetscErrorCode DMPlexTransformView(DMPlexTransform tr, PetscViewer v)
301: {
302:   PetscBool isascii;

304:   PetscFunctionBegin;
306:   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)tr), &v));
308:   PetscCheckSameComm(tr, 1, v, 2);
309:   PetscCall(PetscViewerCheckWritable(v));
310:   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tr, v));
311:   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
312:   if (isascii) PetscCall(DMPlexTransformView_Ascii(tr, v));
313:   PetscTryTypeMethod(tr, view, v);
314:   PetscFunctionReturn(PETSC_SUCCESS);
315: }

317: /*@
318:   DMPlexTransformSetFromOptions - Sets parameters in a transform from values in the options database

320:   Collective

322:   Input Parameter:
323: . tr - the `DMPlexTransform` object to set options for

325:   Options Database Keys:
326: + -dm_plex_transform_type type               - Set the transform type, e.g. refine_regular
327: . -dm_plex_transform_label_match_strata      - Only label points of the same stratum as the producing point
328: . -dm_plex_transform_label_replica_inc inc   - Increment for the label value to be multiplied by the replica number, so that the new label value is oldValue + r * inc
329: . -dm_plex_transform_active name             - Name for active mesh label
330: - -dm_plex_transform_active_values v0,v1,... - Values in the active label

332:   Level: intermediate

334: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformView()`, `DMPlexTransformCreate()`
335: @*/
336: PetscErrorCode DMPlexTransformSetFromOptions(DMPlexTransform tr)
337: {
338:   char        typeName[1024], active[PETSC_MAX_PATH_LEN];
339:   const char *defName = DMPLEXREFINEREGULAR;
340:   PetscBool   flg, match;

342:   PetscFunctionBegin;
344:   PetscObjectOptionsBegin((PetscObject)tr);
345:   PetscCall(PetscOptionsFList("-dm_plex_transform_type", "DMPlexTransform", "DMPlexTransformSetType", DMPlexTransformList, defName, typeName, 1024, &flg));
346:   if (flg) PetscCall(DMPlexTransformSetType(tr, typeName));
347:   else if (!((PetscObject)tr)->type_name) PetscCall(DMPlexTransformSetType(tr, defName));
348:   PetscCall(PetscOptionsBool("-dm_plex_transform_label_match_strata", "Only label points of the same stratum as the producing point", "", tr->labelMatchStrata, &match, &flg));
349:   if (flg) PetscCall(DMPlexTransformSetMatchStrata(tr, match));
350:   PetscCall(PetscOptionsInt("-dm_plex_transform_label_replica_inc", "Increment for the label value to be multiplied by the replica number", "", tr->labelReplicaInc, &tr->labelReplicaInc, NULL));
351:   PetscCall(PetscOptionsString("-dm_plex_transform_active", "Name for active mesh label", "DMPlexTransformSetActive", active, active, sizeof(active), &flg));
352:   if (flg) {
353:     DM       dm;
354:     DMLabel  label;
355:     PetscInt values[16];
356:     PetscInt n = 16;

358:     PetscCall(DMPlexTransformGetDM(tr, &dm));
359:     PetscCall(DMGetLabel(dm, active, &label));
360:     PetscCall(PetscOptionsIntArray("-dm_plex_transform_active_values", "The label values to be active", "DMPlexTransformSetActive", values, &n, &flg));
361:     if (flg && n) {
362:       DMLabel newlabel;

364:       PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Active", &newlabel));
365:       for (PetscInt i = 0; i < n; ++i) {
366:         IS is;

368:         PetscCall(DMLabelGetStratumIS(label, values[i], &is));
369:         PetscCall(DMLabelInsertIS(newlabel, is, values[i]));
370:         PetscCall(ISDestroy(&is));
371:       }
372:       PetscCall(DMPlexTransformSetActive(tr, newlabel));
373:       PetscCall(DMLabelDestroy(&newlabel));
374:     } else {
375:       PetscCall(DMPlexTransformSetActive(tr, label));
376:     }
377:   }
378:   PetscTryTypeMethod(tr, setfromoptions, PetscOptionsObject);
379:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
380:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tr, PetscOptionsObject));
381:   PetscOptionsEnd();
382:   PetscFunctionReturn(PETSC_SUCCESS);
383: }

385: /*@
386:   DMPlexTransformDestroy - Destroys a `DMPlexTransform`

388:   Collective

390:   Input Parameter:
391: . tr - the transform object to destroy

393:   Level: beginner

395: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformView()`, `DMPlexTransformCreate()`
396: @*/
397: PetscErrorCode DMPlexTransformDestroy(DMPlexTransform *tr)
398: {
399:   PetscInt c;

401:   PetscFunctionBegin;
402:   if (!*tr) PetscFunctionReturn(PETSC_SUCCESS);
404:   if (--((PetscObject)*tr)->refct > 0) {
405:     *tr = NULL;
406:     PetscFunctionReturn(PETSC_SUCCESS);
407:   }

409:   PetscTryTypeMethod(*tr, destroy);
410:   PetscCall(DMDestroy(&(*tr)->dm));
411:   PetscCall(DMLabelDestroy(&(*tr)->active));
412:   PetscCall(DMLabelDestroy(&(*tr)->trType));
413:   PetscCall(PetscFree2((*tr)->ctOrderOld, (*tr)->ctOrderInvOld));
414:   PetscCall(PetscFree2((*tr)->ctOrderNew, (*tr)->ctOrderInvNew));
415:   PetscCall(PetscFree2((*tr)->ctStart, (*tr)->ctStartNew));
416:   PetscCall(PetscFree((*tr)->offset));
417:   PetscCall(PetscFree2((*tr)->depthStart, (*tr)->depthEnd));
418:   for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
419:     PetscCall(PetscFEDestroy(&(*tr)->coordFE[c]));
420:     PetscCall(PetscFEGeomDestroy(&(*tr)->refGeom[c]));
421:   }
422:   if ((*tr)->trVerts) {
423:     for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
424:       DMPolytopeType *rct;
425:       PetscInt       *rsize, *rcone, *rornt, Nct, n, r;

427:       if (DMPolytopeTypeGetDim((DMPolytopeType)c) > 0 && c != DM_POLYTOPE_UNKNOWN_CELL && c != DM_POLYTOPE_UNKNOWN_FACE) {
428:         PetscCall(DMPlexTransformCellTransform(*tr, (DMPolytopeType)c, 0, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
429:         for (n = 0; n < Nct; ++n) {
430:           if (rct[n] == DM_POLYTOPE_POINT) continue;
431:           for (r = 0; r < rsize[n]; ++r) PetscCall(PetscFree((*tr)->trSubVerts[c][rct[n]][r]));
432:           PetscCall(PetscFree((*tr)->trSubVerts[c][rct[n]]));
433:         }
434:       }
435:       PetscCall(PetscFree((*tr)->trSubVerts[c]));
436:       PetscCall(PetscFree((*tr)->trVerts[c]));
437:     }
438:   }
439:   PetscCall(PetscFree3((*tr)->trNv, (*tr)->trVerts, (*tr)->trSubVerts));
440:   PetscCall(PetscFree2((*tr)->coordFE, (*tr)->refGeom));
441:   /* We do not destroy (*dm)->data here so that we can reference count backend objects */
442:   PetscCall(PetscHeaderDestroy(tr));
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: static PetscErrorCode DMPlexTransformCreateOffset_Internal(DMPlexTransform tr, PetscInt ctOrderOld[], PetscInt ctStart[], PetscInt **offset)
447: {
448:   DMLabel  trType = tr->trType;
449:   PetscInt c, cN, *off;

451:   PetscFunctionBegin;
452:   if (trType) {
453:     DM              dm;
454:     IS              rtIS;
455:     const PetscInt *reftypes;
456:     PetscInt        Nrt, r;

458:     PetscCall(DMPlexTransformGetDM(tr, &dm));
459:     PetscCall(DMLabelGetNumValues(trType, &Nrt));
460:     PetscCall(DMLabelGetValueIS(trType, &rtIS));
461:     PetscCall(ISGetIndices(rtIS, &reftypes));
462:     PetscCall(PetscCalloc1(Nrt * DM_NUM_POLYTOPES, &off));
463:     for (r = 0; r < Nrt; ++r) {
464:       const PetscInt  rt = reftypes[r];
465:       IS              rtIS;
466:       const PetscInt *points;
467:       DMPolytopeType  ct;
468:       PetscInt        np, p;

470:       PetscCall(DMLabelGetStratumIS(trType, rt, &rtIS));
471:       PetscCall(ISGetLocalSize(rtIS, &np));
472:       PetscCall(ISGetIndices(rtIS, &points));
473:       if (!np) continue;
474:       p = points[0];
475:       PetscCall(ISRestoreIndices(rtIS, &points));
476:       PetscCall(ISDestroy(&rtIS));
477:       PetscCall(DMPlexGetCellType(dm, p, &ct));
478:       for (cN = DM_POLYTOPE_POINT; cN < DM_NUM_POLYTOPES; ++cN) {
479:         const DMPolytopeType ctNew = (DMPolytopeType)cN;
480:         DMPolytopeType      *rct;
481:         PetscInt            *rsize, *cone, *ornt;
482:         PetscInt             Nct, n, s;

484:         if (DMPolytopeTypeGetDim(ct) < 0 || DMPolytopeTypeGetDim(ctNew) < 0) {
485:           off[r * DM_NUM_POLYTOPES + ctNew] = -1;
486:           break;
487:         }
488:         off[r * DM_NUM_POLYTOPES + ctNew] = 0;
489:         for (s = 0; s <= r; ++s) {
490:           const PetscInt st = reftypes[s];
491:           DMPolytopeType sct;
492:           PetscInt       q, qrt;

494:           PetscCall(DMLabelGetStratumIS(trType, st, &rtIS));
495:           PetscCall(ISGetLocalSize(rtIS, &np));
496:           PetscCall(ISGetIndices(rtIS, &points));
497:           if (!np) continue;
498:           q = points[0];
499:           PetscCall(ISRestoreIndices(rtIS, &points));
500:           PetscCall(ISDestroy(&rtIS));
501:           PetscCall(DMPlexGetCellType(dm, q, &sct));
502:           PetscCall(DMPlexTransformCellTransform(tr, sct, q, &qrt, &Nct, &rct, &rsize, &cone, &ornt));
503:           PetscCheck(st == qrt, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Refine type %" PetscInt_FMT " of point %" PetscInt_FMT " does not match predicted type %" PetscInt_FMT, qrt, q, st);
504:           if (st == rt) {
505:             for (n = 0; n < Nct; ++n)
506:               if (rct[n] == ctNew) break;
507:             if (n == Nct) off[r * DM_NUM_POLYTOPES + ctNew] = -1;
508:             break;
509:           }
510:           for (n = 0; n < Nct; ++n) {
511:             if (rct[n] == ctNew) {
512:               PetscInt sn;

514:               PetscCall(DMLabelGetStratumSize(trType, st, &sn));
515:               off[r * DM_NUM_POLYTOPES + ctNew] += sn * rsize[n];
516:             }
517:           }
518:         }
519:       }
520:     }
521:     PetscCall(ISRestoreIndices(rtIS, &reftypes));
522:     PetscCall(ISDestroy(&rtIS));
523:   } else {
524:     PetscCall(PetscCalloc1(DM_NUM_POLYTOPES * DM_NUM_POLYTOPES, &off));
525:     for (c = DM_POLYTOPE_POINT; c < DM_NUM_POLYTOPES; ++c) {
526:       const DMPolytopeType ct = (DMPolytopeType)c;
527:       for (cN = DM_POLYTOPE_POINT; cN < DM_NUM_POLYTOPES; ++cN) {
528:         const DMPolytopeType ctNew = (DMPolytopeType)cN;
529:         DMPolytopeType      *rct;
530:         PetscInt            *rsize, *cone, *ornt;
531:         PetscInt             Nct, n, i;

533:         if (DMPolytopeTypeGetDim(ct) < 0 || ct == DM_POLYTOPE_UNKNOWN_CELL || ct == DM_POLYTOPE_UNKNOWN_FACE || DMPolytopeTypeGetDim(ctNew) < 0 || ctNew == DM_POLYTOPE_UNKNOWN_CELL || ctNew == DM_POLYTOPE_UNKNOWN_FACE) {
534:           off[ct * DM_NUM_POLYTOPES + ctNew] = -1;
535:           continue;
536:         }
537:         off[ct * DM_NUM_POLYTOPES + ctNew] = 0;
538:         for (i = DM_POLYTOPE_POINT; i < DM_NUM_POLYTOPES; ++i) {
539:           const DMPolytopeType ict  = (DMPolytopeType)ctOrderOld[i];
540:           const DMPolytopeType ictn = (DMPolytopeType)ctOrderOld[i + 1];

542:           PetscCall(DMPlexTransformCellTransform(tr, ict, PETSC_DETERMINE, NULL, &Nct, &rct, &rsize, &cone, &ornt));
543:           if (ict == ct) {
544:             for (n = 0; n < Nct; ++n)
545:               if (rct[n] == ctNew) break;
546:             if (n == Nct) off[ct * DM_NUM_POLYTOPES + ctNew] = -1;
547:             break;
548:           }
549:           for (n = 0; n < Nct; ++n)
550:             if (rct[n] == ctNew) off[ct * DM_NUM_POLYTOPES + ctNew] += (ctStart[ictn] - ctStart[ict]) * rsize[n];
551:         }
552:       }
553:     }
554:   }
555:   *offset = off;
556:   PetscFunctionReturn(PETSC_SUCCESS);
557: }

559: /*@
560:   DMPlexTransformSetUp - Create the tables that drive the transform

562:   Input Parameter:
563: . tr - The `DMPlexTransform` object

565:   Level: intermediate

567: .seealso: [](plex_transform_table), [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
568: @*/
569: PetscErrorCode DMPlexTransformSetUp(DMPlexTransform tr)
570: {
571:   DMPolytopeType ctCell;
572:   DM             dm;
573:   PetscInt       pStart, pEnd, p, c, celldim = 0;

575:   PetscFunctionBegin;
577:   if (tr->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
578:   PetscCall(DMPlexTransformGetDM(tr, &dm));
579:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_SetUp, tr, dm, 0, 0));
580:   PetscTryTypeMethod(tr, setup);
581:   PetscCall(DMSetSnapToGeomModel(dm, NULL));
582:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));

584:   if (pEnd > pStart) {
585:     // Ignore cells hanging off of embedded surfaces
586:     PetscInt c = pStart;

588:     ctCell = DM_POLYTOPE_FV_GHOST;
589:     while (DMPolytopeTypeGetDim(ctCell) < 0) PetscCall(DMPlexGetCellType(dm, c++, &ctCell));
590:   } else {
591:     PetscInt dim;

593:     PetscCall(DMGetDimension(dm, &dim));
594:     switch (dim) {
595:     case 0:
596:       ctCell = DM_POLYTOPE_POINT;
597:       break;
598:     case 1:
599:       ctCell = DM_POLYTOPE_SEGMENT;
600:       break;
601:     case 2:
602:       ctCell = DM_POLYTOPE_TRIANGLE;
603:       break;
604:     case 3:
605:       ctCell = DM_POLYTOPE_TETRAHEDRON;
606:       break;
607:     default:
608:       ctCell = DM_POLYTOPE_UNKNOWN;
609:     }
610:   }
611:   PetscCall(DMPlexCreateCellTypeOrder_Internal(dm, DMPolytopeTypeGetDim(ctCell), &tr->ctOrderOld, &tr->ctOrderInvOld));
612:   for (p = pStart; p < pEnd; ++p) {
613:     DMPolytopeType  ct;
614:     DMPolytopeType *rct;
615:     PetscInt       *rsize, *cone, *ornt;
616:     PetscInt        Nct, n;

618:     PetscCall(DMPlexGetCellType(dm, p, &ct));
619:     PetscCheck(ct != DM_POLYTOPE_UNKNOWN && ct != DM_POLYTOPE_UNKNOWN_CELL && ct != DM_POLYTOPE_UNKNOWN_FACE, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No cell type for point %" PetscInt_FMT, p);
620:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &cone, &ornt));
621:     for (n = 0; n < Nct; ++n) celldim = PetscMax(celldim, DMPolytopeTypeGetDim(rct[n]));
622:   }
623:   PetscCall(DMPlexCreateCellTypeOrder_Internal(NULL, celldim, &tr->ctOrderNew, &tr->ctOrderInvNew));
624:   /* Construct sizes and offsets for each cell type */
625:   if (!tr->ctStart) {
626:     PetscInt *ctS, *ctSN, *ctC, *ctCN;

628:     PetscCall(PetscCalloc2(DM_NUM_POLYTOPES + 1, &ctS, DM_NUM_POLYTOPES + 1, &ctSN));
629:     PetscCall(PetscCalloc2(DM_NUM_POLYTOPES + 1, &ctC, DM_NUM_POLYTOPES + 1, &ctCN));
630:     for (p = pStart; p < pEnd; ++p) {
631:       DMPolytopeType  ct;
632:       DMPolytopeType *rct;
633:       PetscInt       *rsize, *cone, *ornt;
634:       PetscInt        Nct, n;

636:       PetscCall(DMPlexGetCellType(dm, p, &ct));
637:       PetscCheck(ct != DM_POLYTOPE_UNKNOWN && ct != DM_POLYTOPE_UNKNOWN_CELL && ct != DM_POLYTOPE_UNKNOWN_FACE, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No cell type for point %" PetscInt_FMT, p);
638:       ++ctC[ct];
639:       PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &cone, &ornt));
640:       for (n = 0; n < Nct; ++n) ctCN[rct[n]] += rsize[n];
641:     }
642:     for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
643:       const PetscInt cto  = tr->ctOrderOld[c];
644:       const PetscInt cton = tr->ctOrderOld[c + 1];
645:       const PetscInt ctn  = tr->ctOrderNew[c];
646:       const PetscInt ctnn = tr->ctOrderNew[c + 1];

648:       ctS[cton]  = ctS[cto] + ctC[cto];
649:       ctSN[ctnn] = ctSN[ctn] + ctCN[ctn];
650:     }
651:     PetscCall(PetscFree2(ctC, ctCN));
652:     tr->ctStart    = ctS;
653:     tr->ctStartNew = ctSN;
654:   }
655:   PetscCall(DMPlexTransformCreateOffset_Internal(tr, tr->ctOrderOld, tr->ctStart, &tr->offset));
656:   // Compute depth information
657:   tr->depth = -1;
658:   for (c = 0; c < DM_NUM_POLYTOPES; ++c)
659:     if (tr->ctStartNew[tr->ctOrderNew[c + 1]] > tr->ctStartNew[tr->ctOrderNew[c]]) tr->depth = PetscMax(tr->depth, DMPolytopeTypeGetDim((DMPolytopeType)tr->ctOrderNew[c]));
660:   PetscCall(PetscMalloc2(tr->depth + 1, &tr->depthStart, tr->depth + 1, &tr->depthEnd));
661:   for (PetscInt d = 0; d <= tr->depth; ++d) {
662:     tr->depthStart[d] = PETSC_INT_MAX;
663:     tr->depthEnd[d]   = -1;
664:   }
665:   for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
666:     const PetscInt dep = DMPolytopeTypeGetDim((DMPolytopeType)tr->ctOrderNew[c]);

668:     if (tr->ctStartNew[tr->ctOrderNew[c + 1]] <= tr->ctStartNew[tr->ctOrderNew[c]]) continue;
669:     tr->depthStart[dep] = PetscMin(tr->depthStart[dep], tr->ctStartNew[tr->ctOrderNew[c]]);
670:     tr->depthEnd[dep]   = PetscMax(tr->depthEnd[dep], tr->ctStartNew[tr->ctOrderNew[c + 1]]);
671:   }
672:   tr->setupcalled = PETSC_TRUE;
673:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_SetUp, tr, dm, 0, 0));
674:   PetscFunctionReturn(PETSC_SUCCESS);
675: }

677: /*@
678:   DMPlexTransformGetDM - Get the base `DM` for the transform

680:   Input Parameter:
681: . tr - The `DMPlexTransform` object

683:   Output Parameter:
684: . dm - The original `DM` which will be transformed

686:   Level: intermediate

688: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformSetDM()`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
689: @*/
690: PetscErrorCode DMPlexTransformGetDM(DMPlexTransform tr, DM *dm)
691: {
692:   PetscFunctionBegin;
694:   PetscAssertPointer(dm, 2);
695:   *dm = tr->dm;
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: /*@
700:   DMPlexTransformSetDM - Set the base `DM` for the transform

702:   Input Parameters:
703: + tr - The `DMPlexTransform` object
704: - dm - The original `DM` which will be transformed

706:   Level: intermediate

708:   Note:
709:   The user does not typically call this, as it is called by `DMPlexTransformApply()`.

711: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformGetDM()`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
712: @*/
713: PetscErrorCode DMPlexTransformSetDM(DMPlexTransform tr, DM dm)
714: {
715:   PetscFunctionBegin;
717:   if (dm) {
719:     PetscCall(PetscObjectReference((PetscObject)dm));
720:   }
721:   PetscCall(DMDestroy(&tr->dm));
722:   tr->dm = dm;
723:   PetscFunctionReturn(PETSC_SUCCESS);
724: }

726: /*@
727:   DMPlexTransformGetActive - Get the `DMLabel` marking the active points for the transform

729:   Input Parameter:
730: . tr - The `DMPlexTransform` object

732:   Output Parameter:
733: . active - The `DMLabel` indicating which points will be transformed

735:   Level: intermediate

737: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformSetActive()`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
738: @*/
739: PetscErrorCode DMPlexTransformGetActive(DMPlexTransform tr, DMLabel *active)
740: {
741:   PetscFunctionBegin;
743:   PetscAssertPointer(active, 2);
744:   *active = tr->active;
745:   PetscFunctionReturn(PETSC_SUCCESS);
746: }

748: /*@
749:   DMPlexTransformSetActive - Set the `DMLabel` marking the active points for the transform

751:   Input Parameters:
752: + tr     - The `DMPlexTransform` object
753: - active - The `DMLabel` indicating which points will be transformed

755:   Level: intermediate

757:   Note:
758:   This only applies to transforms listed in [](plex_transform_table) that operate on a subset of the mesh.

760: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformGetActive()`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
761: @*/
762: PetscErrorCode DMPlexTransformSetActive(DMPlexTransform tr, DMLabel active)
763: {
764:   PetscFunctionBegin;
767:   PetscCall(PetscObjectReference((PetscObject)active));
768:   PetscCall(DMLabelDestroy(&tr->active));
769:   tr->active = active;
770:   PetscFunctionReturn(PETSC_SUCCESS);
771: }

773: /*@
774:   DMPlexTransformGetTransformTypes - Get the `DMLabel` marking the transform type of each point for the transform

776:   Input Parameter:
777: . tr - The `DMPlexTransform` object

779:   Output Parameter:
780: . trType - The `DMLabel` indicating the transform type for each point

782:   Level: intermediate

784: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexSetTransformType()`, `DMPlexTransformGetActive()`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
785: @*/
786: PetscErrorCode DMPlexTransformGetTransformTypes(DMPlexTransform tr, DMLabel *trType)
787: {
788:   PetscFunctionBegin;
790:   PetscAssertPointer(trType, 2);
791:   *trType = tr->trType;
792:   PetscFunctionReturn(PETSC_SUCCESS);
793: }

795: /*@
796:   DMPlexTransformSetTransformTypes - Set the `DMLabel` marking the transform type of each point for the transform

798:   Input Parameters:
799: + tr     - The `DMPlexTransform` object
800: - trType - The original `DM` which will be transformed

802:   Level: intermediate

804: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformGetTransformTypes()`, `DMPlexTransformGetActive())`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
805: @*/
806: PetscErrorCode DMPlexTransformSetTransformTypes(DMPlexTransform tr, DMLabel trType)
807: {
808:   PetscFunctionBegin;
811:   PetscCall(PetscObjectReference((PetscObject)trType));
812:   PetscCall(DMLabelDestroy(&tr->trType));
813:   tr->trType = trType;
814:   PetscFunctionReturn(PETSC_SUCCESS);
815: }

817: static PetscErrorCode DMPlexTransformGetCoordinateFE(DMPlexTransform tr, DMPolytopeType ct, PetscFE *fe)
818: {
819:   PetscFunctionBegin;
820:   if (!tr->coordFE[ct]) {
821:     PetscInt dim, cdim;

823:     dim = DMPolytopeTypeGetDim(ct);
824:     PetscCall(DMGetCoordinateDim(tr->dm, &cdim));
825:     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, cdim, ct, 1, PETSC_DETERMINE, &tr->coordFE[ct]));
826:     {
827:       PetscDualSpace  dsp;
828:       PetscQuadrature quad;
829:       DM              K;
830:       PetscFEGeom    *cg;
831:       PetscScalar    *Xq;
832:       PetscReal      *xq, *wq;
833:       PetscInt        Nq, q;

835:       PetscCall(DMPlexTransformGetCellVertices(tr, ct, &Nq, &Xq));
836:       PetscCall(PetscMalloc1(Nq * cdim, &xq));
837:       for (q = 0; q < Nq * cdim; ++q) xq[q] = PetscRealPart(Xq[q]);
838:       PetscCall(PetscMalloc1(Nq, &wq));
839:       for (q = 0; q < Nq; ++q) wq[q] = 1.0;
840:       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad));
841:       PetscCall(PetscQuadratureSetData(quad, dim, 1, Nq, xq, wq));
842:       PetscCall(PetscFESetQuadrature(tr->coordFE[ct], quad));

844:       PetscCall(PetscFEGetDualSpace(tr->coordFE[ct], &dsp));
845:       PetscCall(PetscDualSpaceGetDM(dsp, &K));
846:       PetscCall(PetscFEGeomCreate(quad, 1, cdim, PETSC_FEGEOM_BASIC, &tr->refGeom[ct]));
847:       cg = tr->refGeom[ct];
848:       PetscCall(DMPlexComputeCellGeometryFEM(K, 0, NULL, cg->v, cg->J, cg->invJ, cg->detJ));
849:       PetscCall(PetscQuadratureDestroy(&quad));
850:     }
851:   }
852:   *fe = tr->coordFE[ct];
853:   PetscFunctionReturn(PETSC_SUCCESS);
854: }

856: PetscErrorCode DMPlexTransformSetDimensions_Internal(DMPlexTransform tr, DM dm, DM tdm)
857: {
858:   PetscInt dim, cdim;

860:   PetscFunctionBegin;
861:   PetscCall(DMGetDimension(dm, &dim));
862:   PetscCall(DMSetDimension(tdm, dim));
863:   PetscCall(DMGetCoordinateDim(dm, &cdim));
864:   PetscCall(DMSetCoordinateDim(tdm, cdim));
865:   PetscFunctionReturn(PETSC_SUCCESS);
866: }

868: /*@
869:   DMPlexTransformSetDimensions - Set the dimensions for the transformed `DM`

871:   Input Parameters:
872: + tr - The `DMPlexTransform` object
873: - dm - The original `DM`

875:   Output Parameter:
876: . trdm - The transformed `DM`

878:   Level: advanced

880: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
881: @*/
882: PetscErrorCode DMPlexTransformSetDimensions(DMPlexTransform tr, DM dm, DM trdm)
883: {
884:   PetscFunctionBegin;
885:   PetscUseTypeMethod(tr, setdimensions, dm, trdm);
886:   PetscFunctionReturn(PETSC_SUCCESS);
887: }

889: PetscErrorCode DMPlexTransformGetChart(DMPlexTransform tr, PetscInt *pStart, PetscInt *pEnd)
890: {
891:   PetscFunctionBegin;
892:   if (pStart) *pStart = 0;
893:   if (pEnd) *pEnd = tr->ctStartNew[tr->ctOrderNew[DM_NUM_POLYTOPES]];
894:   PetscFunctionReturn(PETSC_SUCCESS);
895: }

897: PetscErrorCode DMPlexTransformGetCellType(DMPlexTransform tr, PetscInt cell, DMPolytopeType *celltype)
898: {
899:   PetscInt ctNew;

901:   PetscFunctionBegin;
903:   PetscAssertPointer(celltype, 3);
904:   /* TODO Can do bisection since everything is sorted */
905:   for (ctNew = DM_POLYTOPE_POINT; ctNew < DM_NUM_POLYTOPES; ++ctNew) {
906:     PetscInt ctSN = tr->ctStartNew[ctNew], ctEN = tr->ctStartNew[tr->ctOrderNew[tr->ctOrderInvNew[ctNew] + 1]];

908:     if (cell >= ctSN && cell < ctEN) break;
909:   }
910:   PetscCheck(ctNew < DM_NUM_POLYTOPES, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " cannot be located in the transformed mesh", cell);
911:   *celltype = (DMPolytopeType)ctNew;
912:   PetscFunctionReturn(PETSC_SUCCESS);
913: }

915: PetscErrorCode DMPlexTransformGetCellTypeStratum(DMPlexTransform tr, DMPolytopeType celltype, PetscInt *start, PetscInt *end)
916: {
917:   PetscFunctionBegin;
919:   if (start) *start = tr->ctStartNew[celltype];
920:   if (end) *end = tr->ctStartNew[tr->ctOrderNew[tr->ctOrderInvNew[celltype] + 1]];
921:   PetscFunctionReturn(PETSC_SUCCESS);
922: }

924: PetscErrorCode DMPlexTransformGetDepth(DMPlexTransform tr, PetscInt *depth)
925: {
926:   PetscFunctionBegin;
928:   *depth = tr->depth;
929:   PetscFunctionReturn(PETSC_SUCCESS);
930: }

932: PetscErrorCode DMPlexTransformGetDepthStratum(DMPlexTransform tr, PetscInt depth, PetscInt *start, PetscInt *end)
933: {
934:   PetscFunctionBegin;
936:   if (start) *start = tr->depthStart[depth];
937:   if (end) *end = tr->depthEnd[depth];
938:   PetscFunctionReturn(PETSC_SUCCESS);
939: }

941: /*@
942:   DMPlexTransformGetMatchStrata - Get the flag which determines what points get added to the transformed labels

944:   Not Collective

946:   Input Parameter:
947: . tr - The `DMPlexTransform`

949:   Output Parameter:
950: . match - If `PETSC_TRUE`, only add produced points at the same stratum as the original point to new labels

952:   Level: intermediate

954: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformSetMatchStrata()`, `DMPlexGetPointDepth()`
955: @*/
956: PetscErrorCode DMPlexTransformGetMatchStrata(DMPlexTransform tr, PetscBool *match)
957: {
958:   PetscFunctionBegin;
960:   PetscAssertPointer(match, 2);
961:   *match = tr->labelMatchStrata;
962:   PetscFunctionReturn(PETSC_SUCCESS);
963: }

965: /*@
966:   DMPlexTransformSetMatchStrata - Set the flag which determines what points get added to the transformed labels

968:   Not Collective

970:   Input Parameters:
971: + tr    - The `DMPlexTransform`
972: - match - If `PETSC_TRUE`, only add produced points at the same stratum as the original point to new labels

974:   Level: intermediate

976: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformGetMatchStrata()`, `DMPlexGetPointDepth()`
977: @*/
978: PetscErrorCode DMPlexTransformSetMatchStrata(DMPlexTransform tr, PetscBool match)
979: {
980:   PetscFunctionBegin;
982:   tr->labelMatchStrata = match;
983:   PetscFunctionReturn(PETSC_SUCCESS);
984: }

986: /*@
987:   DMPlexTransformCheck - Verify that the given `DM`, produced by this `DMPlexTransform`, is valid

989:   Input Parameters:
990: + tr - The `DMPlexTransform` object
991: - dm - The `DM` to check

993:   Level: advanced

995: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
996: @*/
997: PetscErrorCode DMPlexTransformCheck(DMPlexTransform tr, DM dm)
998: {
999:   PetscFunctionBegin;
1000:   PetscTryTypeMethod(tr, check, dm);
1001:   PetscFunctionReturn(PETSC_SUCCESS);
1002: }

1004: /*@
1005:   DMPlexTransformGetTargetPoint - Get the number of a point in the transformed mesh based on information from the original mesh.

1007:   Not Collective

1009:   Input Parameters:
1010: + tr    - The `DMPlexTransform`
1011: . ct    - The type of the original point which produces the new point
1012: . ctNew - The type of the new point
1013: . p     - The original point which produces the new point
1014: - r     - The replica number of the new point, meaning it is the rth point of type `ctNew` produced from `p`

1016:   Output Parameter:
1017: . pNew - The new point number

1019:   Level: developer

1021: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformGetSourcePoint()`, `DMPlexTransformCellTransform()`
1022: @*/
1023: PetscErrorCode DMPlexTransformGetTargetPoint(DMPlexTransform tr, DMPolytopeType ct, DMPolytopeType ctNew, PetscInt p, PetscInt r, PetscInt *pNew)
1024: {
1025:   DMPolytopeType *rct;
1026:   PetscInt       *rsize, *cone, *ornt;
1027:   PetscInt        rt, Nct, n, off, rp;
1028:   DMLabel         trType = tr->trType;
1029:   PetscInt        ctS = tr->ctStart[ct], ctE = tr->ctStart[tr->ctOrderOld[tr->ctOrderInvOld[ct] + 1]];
1030:   PetscInt        ctSN = tr->ctStartNew[ctNew], ctEN = tr->ctStartNew[tr->ctOrderNew[tr->ctOrderInvNew[ctNew] + 1]];
1031:   PetscInt        newp = ctSN, cind;

1033:   PetscFunctionBeginHot;
1034:   PetscCheck(p >= ctS && p < ctE, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " is not a %s [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, DMPolytopeTypes[ct], ctS, ctE);
1035:   PetscCall(DMPlexTransformCellTransform(tr, ct, p, &rt, &Nct, &rct, &rsize, &cone, &ornt));
1036:   if (trType) {
1037:     PetscCall(DMLabelGetValueIndex(trType, rt, &cind));
1038:     PetscCall(DMLabelGetStratumPointIndex(trType, rt, p, &rp));
1039:     PetscCheck(rp >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s point %" PetscInt_FMT " does not have refine type %" PetscInt_FMT, DMPolytopeTypes[ct], p, rt);
1040:   } else {
1041:     cind = ct;
1042:     rp   = p - ctS;
1043:   }
1044:   off = tr->offset[cind * DM_NUM_POLYTOPES + ctNew];
1045:   PetscCheck(off >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s (%" PetscInt_FMT ") of point %" PetscInt_FMT " does not produce type %s for transform %s", DMPolytopeTypes[ct], rt, p, DMPolytopeTypes[ctNew], tr->hdr.type_name);
1046:   newp += off;
1047:   for (n = 0; n < Nct; ++n) {
1048:     if (rct[n] == ctNew) {
1049:       if (rsize[n] && r >= rsize[n])
1050:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Replica number %" PetscInt_FMT " for point %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ") for subcell type %s in cell type %s", r, p, rsize[n], DMPolytopeTypes[rct[n]], DMPolytopeTypes[ct]);
1051:       newp += rp * rsize[n] + r;
1052:       break;
1053:     }
1054:   }

1056:   PetscCheck(!(newp < ctSN) && !(newp >= ctEN), PETSC_COMM_SELF, PETSC_ERR_PLIB, "New point %" PetscInt_FMT " is not a %s [%" PetscInt_FMT ", %" PetscInt_FMT ")", newp, DMPolytopeTypes[ctNew], ctSN, ctEN);
1057:   *pNew = newp;
1058:   PetscFunctionReturn(PETSC_SUCCESS);
1059: }

1061: /*@
1062:   DMPlexTransformGetSourcePoint - Get the number of a point in the original mesh based on information from the transformed mesh.

1064:   Not Collective

1066:   Input Parameters:
1067: + tr   - The `DMPlexTransform`
1068: - pNew - The new point number

1070:   Output Parameters:
1071: + ct    - The type of the original point which produces the new point
1072: . ctNew - The type of the new point
1073: . p     - The original point which produces the new point
1074: - r     - The replica number of the new point, meaning it is the rth point of type ctNew produced from p

1076:   Level: developer

1078: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformGetTargetPoint()`, `DMPlexTransformCellTransform()`
1079: @*/
1080: PetscErrorCode DMPlexTransformGetSourcePoint(DMPlexTransform tr, PetscInt pNew, DMPolytopeType *ct, DMPolytopeType *ctNew, PetscInt *p, PetscInt *r)
1081: {
1082:   DMLabel         trType = tr->trType;
1083:   DMPolytopeType *rct, ctN;
1084:   PetscInt       *rsize, *cone, *ornt;
1085:   PetscInt        rt = -1, rtTmp, Nct, n, rp = 0, rO = 0, pO;
1086:   PetscInt        offset = -1, ctS, ctE, ctO = 0, ctTmp, rtS;

1088:   PetscFunctionBegin;
1089:   PetscCall(DMPlexTransformGetCellType(tr, pNew, &ctN));
1090:   if (trType) {
1091:     DM              dm;
1092:     IS              rtIS;
1093:     const PetscInt *reftypes;
1094:     PetscInt        Nrt, r, rtStart;

1096:     PetscCall(DMPlexTransformGetDM(tr, &dm));
1097:     PetscCall(DMLabelGetNumValues(trType, &Nrt));
1098:     PetscCall(DMLabelGetValueIS(trType, &rtIS));
1099:     PetscCall(ISGetIndices(rtIS, &reftypes));
1100:     for (r = 0; r < Nrt; ++r) {
1101:       const PetscInt off = tr->offset[r * DM_NUM_POLYTOPES + ctN];

1103:       if (tr->ctStartNew[ctN] + off > pNew) continue;
1104:       /* Check that any of this refinement type exist */
1105:       /* TODO Actually keep track of the number produced here instead */
1106:       if (off > offset) {
1107:         rt     = reftypes[r];
1108:         offset = off;
1109:       }
1110:     }
1111:     PetscCall(ISRestoreIndices(rtIS, &reftypes));
1112:     PetscCall(ISDestroy(&rtIS));
1113:     PetscCheck(offset >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Source cell type for target point %" PetscInt_FMT " could be not found", pNew);
1114:     /* TODO Map refinement types to cell types */
1115:     PetscCall(DMLabelGetStratumBounds(trType, rt, &rtStart, NULL));
1116:     PetscCheck(rtStart >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Refinement type %" PetscInt_FMT " has no source points", rt);
1117:     for (ctO = 0; ctO < DM_NUM_POLYTOPES; ++ctO) {
1118:       PetscInt ctS = tr->ctStart[ctO], ctE = tr->ctStart[tr->ctOrderOld[tr->ctOrderInvOld[ctO] + 1]];

1120:       if ((rtStart >= ctS) && (rtStart < ctE)) break;
1121:     }
1122:     PetscCheck(ctO != DM_NUM_POLYTOPES, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine a cell type for refinement type %" PetscInt_FMT, rt);
1123:   } else {
1124:     for (ctTmp = 0; ctTmp < DM_NUM_POLYTOPES; ++ctTmp) {
1125:       const PetscInt off = tr->offset[ctTmp * DM_NUM_POLYTOPES + ctN];

1127:       if (tr->ctStartNew[ctN] + off > pNew) continue;
1128:       if (tr->ctStart[tr->ctOrderOld[tr->ctOrderInvOld[ctTmp] + 1]] <= tr->ctStart[ctTmp]) continue;
1129:       /* TODO Actually keep track of the number produced here instead */
1130:       if (off > offset) {
1131:         ctO    = ctTmp;
1132:         offset = off;
1133:       }
1134:     }
1135:     rt = -1;
1136:     PetscCheck(offset >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Source cell type for target point %" PetscInt_FMT " could be not found", pNew);
1137:   }
1138:   ctS = tr->ctStart[ctO];
1139:   ctE = tr->ctStart[tr->ctOrderOld[tr->ctOrderInvOld[ctO] + 1]];
1140:   if (trType) {
1141:     for (rtS = ctS; rtS < ctE; ++rtS) {
1142:       PetscInt val;
1143:       PetscCall(DMLabelGetValue(trType, rtS, &val));
1144:       if (val == rt) break;
1145:     }
1146:     PetscCheck(rtS < ctE, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find point of type %s with refine type %" PetscInt_FMT, DMPolytopeTypes[ctO], rt);
1147:   } else rtS = ctS;
1148:   PetscCall(DMPlexTransformCellTransform(tr, (DMPolytopeType)ctO, rtS, &rtTmp, &Nct, &rct, &rsize, &cone, &ornt));
1149:   PetscCheck(!trType || rt == rtTmp, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " has refine type %" PetscInt_FMT " != %" PetscInt_FMT " refine type which produced point %" PetscInt_FMT, rtS, rtTmp, rt, pNew);
1150:   for (n = 0; n < Nct; ++n) {
1151:     if (rct[n] == ctN) {
1152:       PetscInt tmp = pNew - tr->ctStartNew[ctN] - offset, val, c;

1154:       if (trType) {
1155:         for (c = ctS; c < ctE; ++c) {
1156:           PetscCall(DMLabelGetValue(trType, c, &val));
1157:           if (val == rt) {
1158:             if (tmp < rsize[n]) break;
1159:             tmp -= rsize[n];
1160:           }
1161:         }
1162:         PetscCheck(c < ctE, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent point for target point %" PetscInt_FMT " could be not found", pNew);
1163:         rp = c - ctS;
1164:         rO = tmp;
1165:       } else {
1166:         // This assumes that all points of type ctO transform the same way
1167:         rp = tmp / rsize[n];
1168:         rO = tmp % rsize[n];
1169:       }
1170:       break;
1171:     }
1172:   }
1173:   PetscCheck(n != Nct, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Replica number for target point %" PetscInt_FMT " could be not found", pNew);
1174:   pO = rp + ctS;
1175:   PetscCheck(!(pO < ctS) && !(pO >= ctE), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Source point %" PetscInt_FMT " is not a %s [%" PetscInt_FMT ", %" PetscInt_FMT ")", pO, DMPolytopeTypes[ctO], ctS, ctE);
1176:   if (ct) *ct = (DMPolytopeType)ctO;
1177:   if (ctNew) *ctNew = ctN;
1178:   if (p) *p = pO;
1179:   if (r) *r = rO;
1180:   PetscFunctionReturn(PETSC_SUCCESS);
1181: }

1183: /*@
1184:   DMPlexTransformCellTransform - Describes the transform of a given source cell into a set of other target cells. These produced cells become the new mesh.

1186:   Input Parameters:
1187: + tr     - The `DMPlexTransform` object
1188: . source - The source cell type
1189: - p      - The source point, which can also determine the refine type

1191:   Output Parameters:
1192: + rt     - The refine type for this point
1193: . Nt     - The number of types produced by this point
1194: . target - An array of length `Nt` giving the types produced
1195: . size   - An array of length `Nt` giving the number of cells of each type produced
1196: . cone   - An array of length `Nt`*size[t]*coneSize[t] giving the cell type for each point in the cone of each produced point
1197: - ornt   - An array of length `Nt`*size[t]*coneSize[t] giving the orientation for each point in the cone of each produced point

1199:   Level: advanced

1201:   Notes:
1202:   The cone array gives the cone of each subcell listed by the first three outputs. For each cone point, we
1203:   need the cell type, point identifier, and orientation within the subcell. The orientation is with respect to the canonical
1204:   division (described in these outputs) of the cell in the original mesh. The point identifier is given by
1205: .vb
1206:    the number of cones to be taken, or 0 for the current cell
1207:    the cell cone point number at each level from which it is subdivided
1208:    the replica number r of the subdivision.
1209: .ve
1210:   The orientation is with respect to the canonical cone orientation. For example, the prescription for edge division is
1211: .vb
1212:    Nt     = 2
1213:    target = {DM_POLYTOPE_POINT, DM_POLYTOPE_SEGMENT}
1214:    size   = {1, 2}
1215:    cone   = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 0, 0,  DM_POLYTOPE_POINT, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0}
1216:    ornt   = {                         0,                       0,                        0,                          0}
1217: .ve

1219: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformApply()`, `DMPlexTransformCreate()`
1220: @*/
1221: PetscErrorCode DMPlexTransformCellTransform(DMPlexTransform tr, DMPolytopeType source, PetscInt p, PetscInt *rt, PetscInt *Nt, DMPolytopeType *target[], PetscInt *size[], PetscInt *cone[], PetscInt *ornt[])
1222: {
1223:   PetscFunctionBegin;
1224:   PetscUseTypeMethod(tr, celltransform, source, p, rt, Nt, target, size, cone, ornt);
1225:   PetscFunctionReturn(PETSC_SUCCESS);
1226: }

1228: PetscErrorCode DMPlexTransformGetSubcellOrientationIdentity(DMPlexTransform tr, DMPolytopeType sct, PetscInt sp, PetscInt so, DMPolytopeType tct, PetscInt r, PetscInt o, PetscInt *rnew, PetscInt *onew)
1229: {
1230:   PetscFunctionBegin;
1231:   *rnew = r;
1232:   *onew = DMPolytopeTypeComposeOrientation(tct, o, so);
1233:   PetscFunctionReturn(PETSC_SUCCESS);
1234: }

1236: /* Returns the same thing */
1237: PetscErrorCode DMPlexTransformCellTransformIdentity(DMPlexTransform tr, DMPolytopeType source, PetscInt p, PetscInt *rt, PetscInt *Nt, DMPolytopeType *target[], PetscInt *size[], PetscInt *cone[], PetscInt *ornt[])
1238: {
1239:   static DMPolytopeType vertexT[] = {DM_POLYTOPE_POINT};
1240:   static PetscInt       vertexS[] = {1};
1241:   static PetscInt       vertexC[] = {0};
1242:   static PetscInt       vertexO[] = {0};
1243:   static DMPolytopeType edgeT[]   = {DM_POLYTOPE_SEGMENT};
1244:   static PetscInt       edgeS[]   = {1};
1245:   static PetscInt       edgeC[]   = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0};
1246:   static PetscInt       edgeO[]   = {0, 0};
1247:   static DMPolytopeType tedgeT[]  = {DM_POLYTOPE_POINT_PRISM_TENSOR};
1248:   static PetscInt       tedgeS[]  = {1};
1249:   static PetscInt       tedgeC[]  = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0};
1250:   static PetscInt       tedgeO[]  = {0, 0};
1251:   static DMPolytopeType triT[]    = {DM_POLYTOPE_TRIANGLE};
1252:   static PetscInt       triS[]    = {1};
1253:   static PetscInt       triC[]    = {DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 1, 2, 0};
1254:   static PetscInt       triO[]    = {0, 0, 0};
1255:   static DMPolytopeType quadT[]   = {DM_POLYTOPE_QUADRILATERAL};
1256:   static PetscInt       quadS[]   = {1};
1257:   static PetscInt       quadC[]   = {DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 1, 2, 0, DM_POLYTOPE_SEGMENT, 1, 3, 0};
1258:   static PetscInt       quadO[]   = {0, 0, 0, 0};
1259:   static DMPolytopeType tquadT[]  = {DM_POLYTOPE_SEG_PRISM_TENSOR};
1260:   static PetscInt       tquadS[]  = {1};
1261:   static PetscInt       tquadC[]  = {DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 3, 0};
1262:   static PetscInt       tquadO[]  = {0, 0, 0, 0};
1263:   static DMPolytopeType tetT[]    = {DM_POLYTOPE_TETRAHEDRON};
1264:   static PetscInt       tetS[]    = {1};
1265:   static PetscInt       tetC[]    = {DM_POLYTOPE_TRIANGLE, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_TRIANGLE, 1, 2, 0, DM_POLYTOPE_TRIANGLE, 1, 3, 0};
1266:   static PetscInt       tetO[]    = {0, 0, 0, 0};
1267:   static DMPolytopeType hexT[]    = {DM_POLYTOPE_HEXAHEDRON};
1268:   static PetscInt       hexS[]    = {1};
1269:   static PetscInt       hexC[] = {DM_POLYTOPE_QUADRILATERAL, 1, 0, 0, DM_POLYTOPE_QUADRILATERAL, 1, 1, 0, DM_POLYTOPE_QUADRILATERAL, 1, 2, 0, DM_POLYTOPE_QUADRILATERAL, 1, 3, 0, DM_POLYTOPE_QUADRILATERAL, 1, 4, 0, DM_POLYTOPE_QUADRILATERAL, 1, 5, 0};
1270:   static PetscInt       hexO[] = {0, 0, 0, 0, 0, 0};
1271:   static DMPolytopeType tripT[]   = {DM_POLYTOPE_TRI_PRISM};
1272:   static PetscInt       tripS[]   = {1};
1273:   static PetscInt       tripC[]   = {DM_POLYTOPE_TRIANGLE, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_QUADRILATERAL, 1, 2, 0, DM_POLYTOPE_QUADRILATERAL, 1, 3, 0, DM_POLYTOPE_QUADRILATERAL, 1, 4, 0};
1274:   static PetscInt       tripO[]   = {0, 0, 0, 0, 0};
1275:   static DMPolytopeType ttripT[]  = {DM_POLYTOPE_TRI_PRISM_TENSOR};
1276:   static PetscInt       ttripS[]  = {1};
1277:   static PetscInt       ttripC[]  = {DM_POLYTOPE_TRIANGLE, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 3, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 4, 0};
1278:   static PetscInt       ttripO[]  = {0, 0, 0, 0, 0};
1279:   static DMPolytopeType tquadpT[] = {DM_POLYTOPE_QUAD_PRISM_TENSOR};
1280:   static PetscInt       tquadpS[] = {1};
1281:   static PetscInt       tquadpC[] = {DM_POLYTOPE_QUADRILATERAL,    1, 0, 0, DM_POLYTOPE_QUADRILATERAL,    1, 1, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 2, 0,
1282:                                      DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 3, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 4, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 5, 0};
1283:   static PetscInt       tquadpO[] = {0, 0, 0, 0, 0, 0};
1284:   static DMPolytopeType pyrT[]    = {DM_POLYTOPE_PYRAMID};
1285:   static PetscInt       pyrS[]    = {1};
1286:   static PetscInt       pyrC[]    = {DM_POLYTOPE_QUADRILATERAL, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_TRIANGLE, 1, 2, 0, DM_POLYTOPE_TRIANGLE, 1, 3, 0, DM_POLYTOPE_TRIANGLE, 1, 4, 0};
1287:   static PetscInt       pyrO[]    = {0, 0, 0, 0, 0};

1289:   PetscFunctionBegin;
1290:   if (rt) *rt = 0;
1291:   switch (source) {
1292:   case DM_POLYTOPE_POINT:
1293:     *Nt     = 1;
1294:     *target = vertexT;
1295:     *size   = vertexS;
1296:     *cone   = vertexC;
1297:     *ornt   = vertexO;
1298:     break;
1299:   case DM_POLYTOPE_SEGMENT:
1300:     *Nt     = 1;
1301:     *target = edgeT;
1302:     *size   = edgeS;
1303:     *cone   = edgeC;
1304:     *ornt   = edgeO;
1305:     break;
1306:   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1307:     *Nt     = 1;
1308:     *target = tedgeT;
1309:     *size   = tedgeS;
1310:     *cone   = tedgeC;
1311:     *ornt   = tedgeO;
1312:     break;
1313:   case DM_POLYTOPE_TRIANGLE:
1314:     *Nt     = 1;
1315:     *target = triT;
1316:     *size   = triS;
1317:     *cone   = triC;
1318:     *ornt   = triO;
1319:     break;
1320:   case DM_POLYTOPE_QUADRILATERAL:
1321:     *Nt     = 1;
1322:     *target = quadT;
1323:     *size   = quadS;
1324:     *cone   = quadC;
1325:     *ornt   = quadO;
1326:     break;
1327:   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1328:     *Nt     = 1;
1329:     *target = tquadT;
1330:     *size   = tquadS;
1331:     *cone   = tquadC;
1332:     *ornt   = tquadO;
1333:     break;
1334:   case DM_POLYTOPE_TETRAHEDRON:
1335:     *Nt     = 1;
1336:     *target = tetT;
1337:     *size   = tetS;
1338:     *cone   = tetC;
1339:     *ornt   = tetO;
1340:     break;
1341:   case DM_POLYTOPE_HEXAHEDRON:
1342:     *Nt     = 1;
1343:     *target = hexT;
1344:     *size   = hexS;
1345:     *cone   = hexC;
1346:     *ornt   = hexO;
1347:     break;
1348:   case DM_POLYTOPE_TRI_PRISM:
1349:     *Nt     = 1;
1350:     *target = tripT;
1351:     *size   = tripS;
1352:     *cone   = tripC;
1353:     *ornt   = tripO;
1354:     break;
1355:   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1356:     *Nt     = 1;
1357:     *target = ttripT;
1358:     *size   = ttripS;
1359:     *cone   = ttripC;
1360:     *ornt   = ttripO;
1361:     break;
1362:   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1363:     *Nt     = 1;
1364:     *target = tquadpT;
1365:     *size   = tquadpS;
1366:     *cone   = tquadpC;
1367:     *ornt   = tquadpO;
1368:     break;
1369:   case DM_POLYTOPE_PYRAMID:
1370:     *Nt     = 1;
1371:     *target = pyrT;
1372:     *size   = pyrS;
1373:     *cone   = pyrC;
1374:     *ornt   = pyrO;
1375:     break;
1376:   default:
1377:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No refinement strategy for %s", DMPolytopeTypes[source]);
1378:   }
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

1382: /*@
1383:   DMPlexTransformGetSubcellOrientation - Transform the replica number and orientation for a target point according to the group action for the source point

1385:   Not Collective

1387:   Input Parameters:
1388: + tr  - The `DMPlexTransform`
1389: . sct - The source point cell type, from whom the new cell is being produced
1390: . sp  - The source point
1391: . so  - The orientation of the source point in its enclosing parent
1392: . tct - The target point cell type
1393: . r   - The replica number requested for the produced cell type
1394: - o   - The orientation of the replica

1396:   Output Parameters:
1397: + rnew - The replica number, given the orientation of the parent
1398: - onew - The replica orientation, given the orientation of the parent

1400:   Level: advanced

1402: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformCellTransform()`, `DMPlexTransformApply()`
1403: @*/
1404: PetscErrorCode DMPlexTransformGetSubcellOrientation(DMPlexTransform tr, DMPolytopeType sct, PetscInt sp, PetscInt so, DMPolytopeType tct, PetscInt r, PetscInt o, PetscInt *rnew, PetscInt *onew)
1405: {
1406:   PetscFunctionBeginHot;
1407:   PetscUseTypeMethod(tr, getsubcellorientation, sct, sp, so, tct, r, o, rnew, onew);
1408:   PetscFunctionReturn(PETSC_SUCCESS);
1409: }

1411: static PetscErrorCode DMPlexTransformSetConeSizes(DMPlexTransform tr, DM rdm)
1412: {
1413:   DM       dm;
1414:   PetscInt pStart, pEnd, pNew;

1416:   PetscFunctionBegin;
1417:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1418:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_SetConeSizes, tr, dm, 0, 0));
1419:   /* Must create the celltype label here so that we do not automatically try to compute the types */
1420:   PetscCall(DMCreateLabel(rdm, "celltype"));
1421:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1422:   for (PetscInt p = pStart; p < pEnd; ++p) {
1423:     DMPolytopeType  ct;
1424:     DMPolytopeType *rct;
1425:     PetscInt       *rsize, *rcone, *rornt;
1426:     PetscInt        Nct, n, r;

1428:     PetscCall(DMPlexGetCellType(dm, p, &ct));
1429:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1430:     for (n = 0; n < Nct; ++n) {
1431:       for (r = 0; r < rsize[n]; ++r) {
1432:         PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, r, &pNew));
1433:         PetscCall(DMPlexSetConeSize(rdm, pNew, DMPolytopeTypeGetConeSize(rct[n])));
1434:         PetscCall(DMPlexSetCellType(rdm, pNew, rct[n]));
1435:       }
1436:     }
1437:   }
1438:   /* Let the DM know we have set all the cell types */
1439:   {
1440:     DMLabel  ctLabel;
1441:     DM_Plex *plex = (DM_Plex *)rdm->data;

1443:     PetscCall(DMPlexGetCellTypeLabel(rdm, &ctLabel));
1444:     PetscCall(PetscObjectStateGet((PetscObject)ctLabel, &plex->celltypeState));
1445:   }
1446:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_SetConeSizes, tr, dm, 0, 0));
1447:   PetscFunctionReturn(PETSC_SUCCESS);
1448: }

1450: PetscErrorCode DMPlexTransformGetConeSize(DMPlexTransform tr, PetscInt q, PetscInt *coneSize)
1451: {
1452:   DMPolytopeType ctNew;

1454:   PetscFunctionBegin;
1456:   PetscAssertPointer(coneSize, 3);
1457:   PetscCall(DMPlexTransformGetCellType(tr, q, &ctNew));
1458:   *coneSize = DMPolytopeTypeGetConeSize(ctNew);
1459:   PetscFunctionReturn(PETSC_SUCCESS);
1460: }

1462: /* The orientation o is for the interior of the cell p */
1463: static PetscErrorCode DMPlexTransformGetCone_Internal(DMPlexTransform tr, PetscInt p, PetscInt o, DMPolytopeType ct, DMPolytopeType ctNew, const PetscInt rcone[], PetscInt *coneoff, const PetscInt rornt[], PetscInt *orntoff, PetscInt coneNew[], PetscInt orntNew[])
1464: {
1465:   DM              dm;
1466:   const PetscInt  csizeNew = DMPolytopeTypeGetConeSize(ctNew);
1467:   const PetscInt *cone;
1468:   DMPolytopeType *newft = NULL;
1469:   PetscInt        c, coff = *coneoff, ooff = *orntoff;
1470:   PetscInt        dim, cr = 0, co = 0, nr, no;

1472:   PetscFunctionBegin;
1473:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1474:   PetscCall(DMPlexGetOrientedCone(dm, p, &cone, NULL));
1475:   // Check if we have to permute this cell
1476:   PetscCall(DMGetDimension(dm, &dim));
1477:   if (DMPolytopeTypeGetDim(ctNew) == dim && DMPolytopeTypeGetDim(ct) == dim - 1) {
1478:     PetscCall(DMPlexTransformGetSubcellOrientation(tr, ct, p, o, ctNew, cr, co, &nr, &no));
1479:     if (cr != nr || co != no) PetscCall(DMGetWorkArray(dm, csizeNew, MPIU_INT, &newft));
1480:   }
1481:   for (c = 0; c < csizeNew; ++c) {
1482:     PetscInt             ppp   = -1;                            /* Parent Parent point: Parent of point pp */
1483:     PetscInt             pp    = p;                             /* Parent point: Point in the original mesh producing new cone point */
1484:     PetscInt             po    = 0;                             /* Orientation of parent point pp in parent parent point ppp */
1485:     DMPolytopeType       pct   = ct;                            /* Parent type: Cell type for parent of new cone point */
1486:     const PetscInt      *pcone = cone;                          /* Parent cone: Cone of parent point pp */
1487:     PetscInt             pr    = -1;                            /* Replica number of pp that produces new cone point  */
1488:     const DMPolytopeType ft    = (DMPolytopeType)rcone[coff++]; /* Cell type for new cone point of pNew */
1489:     const PetscInt       fn    = rcone[coff++];                 /* Number of cones of p that need to be taken when producing new cone point */
1490:     PetscInt             fo    = rornt[ooff++];                 /* Orientation of new cone point in pNew */
1491:     PetscInt             lc;

1493:     /* Get the type (pct) and point number (pp) of the parent point in the original mesh which produces this cone point */
1494:     for (lc = 0; lc < fn; ++lc) {
1495:       const PetscInt *parr = DMPolytopeTypeGetArrangement(pct, po);
1496:       const PetscInt  acp  = rcone[coff++];
1497:       const PetscInt  pcp  = parr[acp * 2];
1498:       const PetscInt  pco  = parr[acp * 2 + 1];
1499:       const PetscInt *ppornt;

1501:       ppp = pp;
1502:       pp  = pcone[pcp];
1503:       PetscCall(DMPlexGetCellType(dm, pp, &pct));
1504:       // Restore the parent cone from the last iterate
1505:       if (lc) PetscCall(DMPlexRestoreOrientedCone(dm, ppp, &pcone, NULL));
1506:       PetscCall(DMPlexGetOrientedCone(dm, pp, &pcone, NULL));
1507:       PetscCall(DMPlexGetOrientedCone(dm, ppp, NULL, &ppornt));
1508:       po = DMPolytopeTypeComposeOrientation(pct, ppornt[pcp], pco);
1509:       PetscCall(DMPlexRestoreOrientedCone(dm, ppp, NULL, &ppornt));
1510:     }
1511:     if (lc) PetscCall(DMPlexRestoreOrientedCone(dm, pp, &pcone, NULL));
1512:     pr = rcone[coff++];
1513:     /* Orientation po of pp maps (pr, fo) -> (pr', fo') */
1514:     PetscCall(DMPlexTransformGetSubcellOrientation(tr, pct, pp, fn ? po : o, ft, pr, fo, &pr, &fo));
1515:     PetscCall(DMPlexTransformGetTargetPoint(tr, pct, ft, pp, pr, &coneNew[c]));
1516:     orntNew[c] = fo;
1517:     if (newft) newft[c] = ft;
1518:   }
1519:   PetscCall(DMPlexRestoreOrientedCone(dm, p, &cone, NULL));
1520:   if (newft) {
1521:     const PetscInt *arr;
1522:     PetscInt       *newcone, *newornt;

1524:     arr = DMPolytopeTypeGetArrangement(ctNew, no);
1525:     PetscCall(DMGetWorkArray(dm, csizeNew, MPIU_INT, &newcone));
1526:     PetscCall(DMGetWorkArray(dm, csizeNew, MPIU_INT, &newornt));
1527:     for (PetscInt c = 0; c < csizeNew; ++c) {
1528:       DMPolytopeType ft = newft[c];
1529:       PetscInt       nO;

1531:       nO         = DMPolytopeTypeGetNumArrangements(ft) / 2;
1532:       newcone[c] = coneNew[arr[c * 2 + 0]];
1533:       newornt[c] = DMPolytopeTypeComposeOrientation(ft, arr[c * 2 + 1], orntNew[arr[c * 2 + 0]]);
1534:       PetscCheck(!newornt[c] || !(newornt[c] >= nO || newornt[c] < -nO), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, newornt[c], -nO, nO, DMPolytopeTypes[ft], coneNew[c]);
1535:     }
1536:     for (PetscInt c = 0; c < csizeNew; ++c) {
1537:       coneNew[c] = newcone[c];
1538:       orntNew[c] = newornt[c];
1539:     }
1540:     PetscCall(DMRestoreWorkArray(dm, csizeNew, MPIU_INT, &newcone));
1541:     PetscCall(DMRestoreWorkArray(dm, csizeNew, MPIU_INT, &newornt));
1542:     PetscCall(DMRestoreWorkArray(dm, csizeNew, MPIU_INT, &newft));
1543:   }
1544:   *coneoff = coff;
1545:   *orntoff = ooff;
1546:   PetscFunctionReturn(PETSC_SUCCESS);
1547: }

1549: static PetscErrorCode DMPlexTransformSetCones(DMPlexTransform tr, DM rdm)
1550: {
1551:   DM             dm;
1552:   DMPolytopeType ct;
1553:   PetscInt      *coneNew, *orntNew;
1554:   PetscInt       maxConeSize = 0, pStart, pEnd, p, pNew;

1556:   PetscFunctionBegin;
1557:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1558:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_SetCones, tr, dm, 0, 0));
1559:   for (p = 0; p < DM_NUM_POLYTOPES; ++p) maxConeSize = PetscMax(maxConeSize, DMPolytopeTypeGetConeSize((DMPolytopeType)p));
1560:   PetscCall(DMGetWorkArray(rdm, maxConeSize, MPIU_INT, &coneNew));
1561:   PetscCall(DMGetWorkArray(rdm, maxConeSize, MPIU_INT, &orntNew));
1562:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1563:   for (p = pStart; p < pEnd; ++p) {
1564:     PetscInt        coff, ooff;
1565:     DMPolytopeType *rct;
1566:     PetscInt       *rsize, *rcone, *rornt;
1567:     PetscInt        Nct, n, r;

1569:     PetscCall(DMPlexGetCellType(dm, p, &ct));
1570:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1571:     for (n = 0, coff = 0, ooff = 0; n < Nct; ++n) {
1572:       const DMPolytopeType ctNew = rct[n];

1574:       for (r = 0; r < rsize[n]; ++r) {
1575:         PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, r, &pNew));
1576:         PetscCall(DMPlexTransformGetCone_Internal(tr, p, 0, ct, ctNew, rcone, &coff, rornt, &ooff, coneNew, orntNew));
1577:         PetscCall(DMPlexSetCone(rdm, pNew, coneNew));
1578:         PetscCall(DMPlexSetConeOrientation(rdm, pNew, orntNew));
1579:       }
1580:     }
1581:   }
1582:   PetscCall(DMRestoreWorkArray(rdm, maxConeSize, MPIU_INT, &coneNew));
1583:   PetscCall(DMRestoreWorkArray(rdm, maxConeSize, MPIU_INT, &orntNew));
1584:   PetscCall(DMViewFromOptions(rdm, NULL, "-rdm_view"));
1585:   PetscCall(DMPlexSymmetrize(rdm));
1586:   PetscCall(DMPlexStratify(rdm));
1587:   PetscTryTypeMethod(tr, ordersupports, dm, rdm);
1588:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_SetCones, tr, dm, 0, 0));
1589:   PetscFunctionReturn(PETSC_SUCCESS);
1590: }

1592: PetscErrorCode DMPlexTransformGetConeOriented(DMPlexTransform tr, PetscInt q, PetscInt po, const PetscInt *cone[], const PetscInt *ornt[])
1593: {
1594:   DM              dm;
1595:   DMPolytopeType  ct, qct;
1596:   DMPolytopeType *rct;
1597:   PetscInt       *rsize, *rcone, *rornt, *qcone, *qornt;
1598:   PetscInt        maxConeSize = 0, Nct, p, r, n, nr, coff = 0, ooff = 0;

1600:   PetscFunctionBegin;
1602:   PetscAssertPointer(cone, 4);
1603:   PetscAssertPointer(ornt, 5);
1604:   for (p = 0; p < DM_NUM_POLYTOPES; ++p) maxConeSize = PetscMax(maxConeSize, DMPolytopeTypeGetConeSize((DMPolytopeType)p));
1605:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1606:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &qcone));
1607:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &qornt));
1608:   PetscCall(DMPlexTransformGetSourcePoint(tr, q, &ct, &qct, &p, &r));
1609:   PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1610:   for (n = 0; n < Nct; ++n) {
1611:     const DMPolytopeType ctNew    = rct[n];
1612:     const PetscInt       csizeNew = DMPolytopeTypeGetConeSize(ctNew);
1613:     PetscInt             Nr       = rsize[n], fn, c;

1615:     if (ctNew == qct) Nr = r;
1616:     for (nr = 0; nr < Nr; ++nr) {
1617:       for (c = 0; c < csizeNew; ++c) {
1618:         ++coff;             /* Cell type of new cone point */
1619:         fn = rcone[coff++]; /* Number of cones of p that need to be taken when producing new cone point */
1620:         coff += fn;
1621:         ++coff; /* Replica number of new cone point */
1622:         ++ooff; /* Orientation of new cone point */
1623:       }
1624:     }
1625:     if (ctNew == qct) break;
1626:   }
1627:   PetscCall(DMPlexTransformGetCone_Internal(tr, p, po, ct, qct, rcone, &coff, rornt, &ooff, qcone, qornt));
1628:   *cone = qcone;
1629:   *ornt = qornt;
1630:   PetscFunctionReturn(PETSC_SUCCESS);
1631: }

1633: PetscErrorCode DMPlexTransformGetCone(DMPlexTransform tr, PetscInt q, const PetscInt *cone[], const PetscInt *ornt[])
1634: {
1635:   DM              dm;
1636:   DMPolytopeType  ct, qct;
1637:   DMPolytopeType *rct;
1638:   PetscInt       *rsize, *rcone, *rornt, *qcone, *qornt;
1639:   PetscInt        maxConeSize = 0, Nct, p, r, n, nr, coff = 0, ooff = 0;

1641:   PetscFunctionBegin;
1643:   if (cone) PetscAssertPointer(cone, 3);
1644:   if (ornt) PetscAssertPointer(ornt, 4);
1645:   for (p = 0; p < DM_NUM_POLYTOPES; ++p) maxConeSize = PetscMax(maxConeSize, DMPolytopeTypeGetConeSize((DMPolytopeType)p));
1646:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1647:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &qcone));
1648:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &qornt));
1649:   PetscCall(DMPlexTransformGetSourcePoint(tr, q, &ct, &qct, &p, &r));
1650:   PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1651:   for (n = 0; n < Nct; ++n) {
1652:     const DMPolytopeType ctNew    = rct[n];
1653:     const PetscInt       csizeNew = DMPolytopeTypeGetConeSize(ctNew);
1654:     PetscInt             Nr       = rsize[n], fn, c;

1656:     if (ctNew == qct) Nr = r;
1657:     for (nr = 0; nr < Nr; ++nr) {
1658:       for (c = 0; c < csizeNew; ++c) {
1659:         ++coff;             /* Cell type of new cone point */
1660:         fn = rcone[coff++]; /* Number of cones of p that need to be taken when producing new cone point */
1661:         coff += fn;
1662:         ++coff; /* Replica number of new cone point */
1663:         ++ooff; /* Orientation of new cone point */
1664:       }
1665:     }
1666:     if (ctNew == qct) break;
1667:   }
1668:   PetscCall(DMPlexTransformGetCone_Internal(tr, p, 0, ct, qct, rcone, &coff, rornt, &ooff, qcone, qornt));
1669:   if (cone) *cone = qcone;
1670:   else PetscCall(DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &qcone));
1671:   if (ornt) *ornt = qornt;
1672:   else PetscCall(DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &qornt));
1673:   PetscFunctionReturn(PETSC_SUCCESS);
1674: }

1676: PetscErrorCode DMPlexTransformRestoreCone(DMPlexTransform tr, PetscInt q, const PetscInt *cone[], const PetscInt *ornt[])
1677: {
1678:   DM dm;

1680:   PetscFunctionBegin;
1682:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1683:   if (cone) PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, cone));
1684:   if (ornt) PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, ornt));
1685:   PetscFunctionReturn(PETSC_SUCCESS);
1686: }

1688: static PetscErrorCode DMPlexTransformCreateCellVertices_Internal(DMPlexTransform tr)
1689: {
1690:   PetscInt ict;

1692:   PetscFunctionBegin;
1693:   PetscCall(PetscCalloc3(DM_NUM_POLYTOPES, &tr->trNv, DM_NUM_POLYTOPES, &tr->trVerts, DM_NUM_POLYTOPES, &tr->trSubVerts));
1694:   for (ict = DM_POLYTOPE_POINT; ict < DM_NUM_POLYTOPES; ++ict) {
1695:     const DMPolytopeType ct = (DMPolytopeType)ict;
1696:     DMPlexTransform      reftr;
1697:     DM                   refdm, trdm;
1698:     Vec                  coordinates;
1699:     const PetscScalar   *coords;
1700:     DMPolytopeType      *rct;
1701:     PetscInt            *rsize, *rcone, *rornt;
1702:     PetscInt             Nct, n, r, pNew = 0;
1703:     PetscInt             trdim, vStart, vEnd, Nc;
1704:     const PetscInt       debug = 0;
1705:     const char          *typeName;

1707:     /* Since points are 0-dimensional, coordinates make no sense */
1708:     if (DMPolytopeTypeGetDim(ct) <= 0 || ct == DM_POLYTOPE_UNKNOWN_CELL || ct == DM_POLYTOPE_UNKNOWN_FACE) continue;
1709:     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &refdm));
1710:     PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &reftr));
1711:     PetscCall(DMPlexTransformSetDM(reftr, refdm));
1712:     PetscCall(DMPlexTransformGetType(tr, &typeName));
1713:     PetscCall(DMPlexTransformSetType(reftr, typeName));
1714:     PetscCall(DMPlexTransformSetUp(reftr));
1715:     PetscCall(DMPlexTransformApply(reftr, refdm, &trdm));

1717:     PetscCall(DMGetDimension(trdm, &trdim));
1718:     PetscCall(DMPlexGetDepthStratum(trdm, 0, &vStart, &vEnd));
1719:     tr->trNv[ct] = vEnd - vStart;
1720:     PetscCall(DMGetCoordinatesLocal(trdm, &coordinates));
1721:     PetscCall(VecGetLocalSize(coordinates, &Nc));
1722:     PetscCheck(tr->trNv[ct] * trdim == Nc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell type %s, transformed coordinate size %" PetscInt_FMT " != %" PetscInt_FMT " size of coordinate storage", DMPolytopeTypes[ct], tr->trNv[ct] * trdim, Nc);
1723:     PetscCall(PetscCalloc1(Nc, &tr->trVerts[ct]));
1724:     PetscCall(VecGetArrayRead(coordinates, &coords));
1725:     PetscCall(PetscArraycpy(tr->trVerts[ct], coords, Nc));
1726:     PetscCall(VecRestoreArrayRead(coordinates, &coords));

1728:     PetscCall(PetscCalloc1(DM_NUM_POLYTOPES, &tr->trSubVerts[ct]));
1729:     PetscCall(DMPlexTransformCellTransform(reftr, ct, 0, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1730:     for (n = 0; n < Nct; ++n) {
1731:       /* Since points are 0-dimensional, coordinates make no sense */
1732:       if (rct[n] == DM_POLYTOPE_POINT) continue;
1733:       PetscCall(PetscCalloc1(rsize[n], &tr->trSubVerts[ct][rct[n]]));
1734:       for (r = 0; r < rsize[n]; ++r) {
1735:         PetscInt *closure = NULL;
1736:         PetscInt  clSize, cl, Nv = 0;

1738:         PetscCall(PetscCalloc1(DMPolytopeTypeGetNumVertices(rct[n]), &tr->trSubVerts[ct][rct[n]][r]));
1739:         PetscCall(DMPlexTransformGetTargetPoint(reftr, ct, rct[n], 0, r, &pNew));
1740:         PetscCall(DMPlexGetTransitiveClosure(trdm, pNew, PETSC_TRUE, &clSize, &closure));
1741:         for (cl = 0; cl < clSize * 2; cl += 2) {
1742:           const PetscInt sv = closure[cl];

1744:           if ((sv >= vStart) && (sv < vEnd)) tr->trSubVerts[ct][rct[n]][r][Nv++] = sv - vStart;
1745:         }
1746:         PetscCall(DMPlexRestoreTransitiveClosure(trdm, pNew, PETSC_TRUE, &clSize, &closure));
1747:         PetscCheck(Nv == DMPolytopeTypeGetNumVertices(rct[n]), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of vertices %" PetscInt_FMT " != %" PetscInt_FMT " for %s subcell %" PetscInt_FMT " from cell %s", Nv, DMPolytopeTypeGetNumVertices(rct[n]), DMPolytopeTypes[rct[n]], r, DMPolytopeTypes[ct]);
1748:       }
1749:     }
1750:     if (debug) {
1751:       DMPolytopeType *rct;
1752:       PetscInt       *rsize, *rcone, *rornt;
1753:       PetscInt        v, dE = trdim, d, off = 0;

1755:       PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s: %" PetscInt_FMT " vertices\n", DMPolytopeTypes[ct], tr->trNv[ct]));
1756:       for (v = 0; v < tr->trNv[ct]; ++v) {
1757:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  "));
1758:         for (d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g ", (double)PetscRealPart(tr->trVerts[ct][off++])));
1759:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1760:       }

1762:       PetscCall(DMPlexTransformCellTransform(reftr, ct, 0, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1763:       for (n = 0; n < Nct; ++n) {
1764:         if (rct[n] == DM_POLYTOPE_POINT) continue;
1765:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s: %s subvertices %" PetscInt_FMT "\n", DMPolytopeTypes[ct], DMPolytopeTypes[rct[n]], tr->trNv[ct]));
1766:         for (r = 0; r < rsize[n]; ++r) {
1767:           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  "));
1768:           for (v = 0; v < DMPolytopeTypeGetNumVertices(rct[n]); ++v) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " ", tr->trSubVerts[ct][rct[n]][r][v]));
1769:           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1770:         }
1771:       }
1772:     }
1773:     PetscCall(DMDestroy(&refdm));
1774:     PetscCall(DMDestroy(&trdm));
1775:     PetscCall(DMPlexTransformDestroy(&reftr));
1776:   }
1777:   PetscFunctionReturn(PETSC_SUCCESS);
1778: }

1780: /*@C
1781:   DMPlexTransformGetCellVertices - Get the set of transformed vertices lying in the closure of a reference cell of given type

1783:   Input Parameters:
1784: + tr - The `DMPlexTransform` object
1785: - ct - The cell type

1787:   Output Parameters:
1788: + Nv      - The number of transformed vertices in the closure of the reference cell of given type
1789: - trVerts - The coordinates of these vertices in the reference cell

1791:   Level: developer

1793: .seealso: `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformGetSubcellVertices()`
1794: @*/
1795: PetscErrorCode DMPlexTransformGetCellVertices(DMPlexTransform tr, DMPolytopeType ct, PetscInt *Nv, PetscScalar *trVerts[])
1796: {
1797:   PetscFunctionBegin;
1798:   if (!tr->trNv) PetscCall(DMPlexTransformCreateCellVertices_Internal(tr));
1799:   if (Nv) *Nv = tr->trNv[ct];
1800:   if (trVerts) *trVerts = tr->trVerts[ct];
1801:   PetscFunctionReturn(PETSC_SUCCESS);
1802: }

1804: /*@C
1805:   DMPlexTransformGetSubcellVertices - Get the set of transformed vertices defining a subcell in the reference cell of given type

1807:   Input Parameters:
1808: + tr  - The `DMPlexTransform` object
1809: . ct  - The cell type
1810: . rct - The subcell type
1811: - r   - The subcell index

1813:   Output Parameter:
1814: . subVerts - The indices of these vertices in the set of vertices returned by `DMPlexTransformGetCellVertices()`

1816:   Level: developer

1818: .seealso: `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformGetCellVertices()`
1819: @*/
1820: PetscErrorCode DMPlexTransformGetSubcellVertices(DMPlexTransform tr, DMPolytopeType ct, DMPolytopeType rct, PetscInt r, PetscInt *subVerts[])
1821: {
1822:   PetscFunctionBegin;
1823:   if (!tr->trNv) PetscCall(DMPlexTransformCreateCellVertices_Internal(tr));
1824:   PetscCheck(tr->trSubVerts[ct][rct], PetscObjectComm((PetscObject)tr), PETSC_ERR_ARG_WRONG, "Cell type %s does not produce %s", DMPolytopeTypes[ct], DMPolytopeTypes[rct]);
1825:   if (subVerts) *subVerts = tr->trSubVerts[ct][rct][r];
1826:   PetscFunctionReturn(PETSC_SUCCESS);
1827: }

1829: /* Computes new vertex as the barycenter, or centroid */
1830: PetscErrorCode DMPlexTransformMapCoordinatesBarycenter_Internal(DMPlexTransform tr, DMPolytopeType pct, DMPolytopeType ct, PetscInt p, PetscInt r, PetscInt Nv, PetscInt dE, const PetscScalar in[], PetscScalar out[])
1831: {
1832:   PetscInt v, d;

1834:   PetscFunctionBeginHot;
1835:   PetscCheck(ct == DM_POLYTOPE_POINT, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not for refined point type %s", DMPolytopeTypes[ct]);
1836:   for (d = 0; d < dE; ++d) out[d] = 0.0;
1837:   for (v = 0; v < Nv; ++v)
1838:     for (d = 0; d < dE; ++d) out[d] += in[v * dE + d];
1839:   for (d = 0; d < dE; ++d) out[d] /= Nv;
1840:   PetscFunctionReturn(PETSC_SUCCESS);
1841: }

1843: /*@
1844:   DMPlexTransformMapCoordinates - Calculate new coordinates for produced points

1846:   Not collective

1848:   Input Parameters:
1849: + tr  - The `DMPlexTransform`
1850: . pct - The cell type of the parent, from whom the new cell is being produced
1851: . ct  - The type being produced
1852: . p   - The original point
1853: . r   - The replica number requested for the produced cell type
1854: . Nv  - Number of vertices in the closure of the parent cell
1855: . dE  - Spatial dimension
1856: - in  - array of size Nv*dE, holding coordinates of the vertices in the closure of the parent cell

1858:   Output Parameter:
1859: . out - The coordinates of the new vertices

1861:   Level: intermediate

1863: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexTransformApply()`
1864: @*/
1865: PetscErrorCode DMPlexTransformMapCoordinates(DMPlexTransform tr, DMPolytopeType pct, DMPolytopeType ct, PetscInt p, PetscInt r, PetscInt Nv, PetscInt dE, const PetscScalar in[], PetscScalar out[])
1866: {
1867:   PetscFunctionBeginHot;
1868:   if (Nv) PetscUseTypeMethod(tr, mapcoordinates, pct, ct, p, r, Nv, dE, in, out);
1869:   PetscFunctionReturn(PETSC_SUCCESS);
1870: }

1872: /*
1873:   DMPlexTransformLabelProducedPoint_Private - Label a produced point based on its parent label

1875:   Not Collective

1877:   Input Parameters:
1878: + tr    - The `DMPlexTransform`
1879: . label - The label in the transformed mesh
1880: . pp    - The parent point in the original mesh
1881: . pct   - The cell type of the parent point
1882: . p     - The point in the transformed mesh
1883: . ct    - The cell type of the point
1884: . r     - The replica number of the point
1885: - val   - The label value of the parent point

1887:   Level: developer

1889: .seealso: `DMPlexTransformCreateLabels()`, `RefineLabel_Internal()`
1890: */
1891: static PetscErrorCode DMPlexTransformLabelProducedPoint_Private(DMPlexTransform tr, DMLabel label, PetscInt pp, DMPolytopeType pct, PetscInt p, DMPolytopeType ct, PetscInt r, PetscInt val)
1892: {
1893:   PetscFunctionBeginHot;
1894:   if (tr->labelMatchStrata && pct != ct) PetscFunctionReturn(PETSC_SUCCESS);
1895:   PetscCall(DMLabelSetValue(label, p, val + tr->labelReplicaInc * r));
1896:   PetscFunctionReturn(PETSC_SUCCESS);
1897: }

1899: static PetscErrorCode RefineLabel_Internal(DMPlexTransform tr, DMLabel label, DMLabel labelNew)
1900: {
1901:   DM              dm;
1902:   IS              valueIS;
1903:   const PetscInt *values;
1904:   PetscInt        defVal, Nv, val;

1906:   PetscFunctionBegin;
1907:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1908:   PetscCall(DMLabelGetDefaultValue(label, &defVal));
1909:   PetscCall(DMLabelSetDefaultValue(labelNew, defVal));
1910:   PetscCall(DMLabelGetValueIS(label, &valueIS));
1911:   PetscCall(ISGetLocalSize(valueIS, &Nv));
1912:   PetscCall(ISGetIndices(valueIS, &values));
1913:   for (val = 0; val < Nv; ++val) {
1914:     IS              pointIS;
1915:     const PetscInt *points;
1916:     PetscInt        numPoints, p;

1918:     /* Ensure refined label is created with same number of strata as
1919:      * original (even if no entries here). */
1920:     PetscCall(DMLabelAddStratum(labelNew, values[val]));
1921:     PetscCall(DMLabelGetStratumIS(label, values[val], &pointIS));
1922:     PetscCall(ISGetLocalSize(pointIS, &numPoints));
1923:     PetscCall(ISGetIndices(pointIS, &points));
1924:     for (p = 0; p < numPoints; ++p) {
1925:       const PetscInt  point = points[p];
1926:       DMPolytopeType  ct;
1927:       DMPolytopeType *rct;
1928:       PetscInt       *rsize, *rcone, *rornt;
1929:       PetscInt        Nct, n, r, pNew = 0;

1931:       PetscCall(DMPlexGetCellType(dm, point, &ct));
1932:       PetscCall(DMPlexTransformCellTransform(tr, ct, point, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
1933:       for (n = 0; n < Nct; ++n) {
1934:         for (r = 0; r < rsize[n]; ++r) {
1935:           PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], point, r, &pNew));
1936:           PetscCall(DMPlexTransformLabelProducedPoint_Private(tr, labelNew, point, ct, pNew, rct[n], r, values[val]));
1937:         }
1938:       }
1939:     }
1940:     PetscCall(ISRestoreIndices(pointIS, &points));
1941:     PetscCall(ISDestroy(&pointIS));
1942:   }
1943:   PetscCall(ISRestoreIndices(valueIS, &values));
1944:   PetscCall(ISDestroy(&valueIS));
1945:   PetscFunctionReturn(PETSC_SUCCESS);
1946: }

1948: static PetscErrorCode DMPlexTransformCreateLabels(DMPlexTransform tr, DM rdm)
1949: {
1950:   DM       dm;
1951:   PetscInt numLabels, l;

1953:   PetscFunctionBegin;
1954:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1955:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_CreateLabels, tr, dm, 0, 0));
1956:   PetscCall(DMGetNumLabels(dm, &numLabels));
1957:   for (l = 0; l < numLabels; ++l) {
1958:     DMLabel     label, labelNew;
1959:     const char *lname;
1960:     PetscBool   isDepth, isCellType;

1962:     PetscCall(DMGetLabelName(dm, l, &lname));
1963:     PetscCall(PetscStrcmp(lname, "depth", &isDepth));
1964:     if (isDepth) continue;
1965:     PetscCall(PetscStrcmp(lname, "celltype", &isCellType));
1966:     if (isCellType) continue;
1967:     PetscCall(DMCreateLabel(rdm, lname));
1968:     PetscCall(DMGetLabel(dm, lname, &label));
1969:     PetscCall(DMGetLabel(rdm, lname, &labelNew));
1970:     PetscCall(RefineLabel_Internal(tr, label, labelNew));
1971:   }
1972:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_CreateLabels, tr, dm, 0, 0));
1973:   PetscFunctionReturn(PETSC_SUCCESS);
1974: }

1976: /* This refines the labels which define regions for fields and DSes since they are not in the list of labels for the DM */
1977: PetscErrorCode DMPlexTransformCreateDiscLabels(DMPlexTransform tr, DM rdm)
1978: {
1979:   DM       dm;
1980:   PetscInt Nf, f, Nds, s;

1982:   PetscFunctionBegin;
1983:   PetscCall(DMPlexTransformGetDM(tr, &dm));
1984:   PetscCall(DMGetNumFields(dm, &Nf));
1985:   for (f = 0; f < Nf; ++f) {
1986:     DMLabel     label, labelNew;
1987:     PetscObject obj;
1988:     const char *lname;

1990:     PetscCall(DMGetField(rdm, f, &label, &obj));
1991:     if (!label) continue;
1992:     PetscCall(PetscObjectGetName((PetscObject)label, &lname));
1993:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, lname, &labelNew));
1994:     PetscCall(RefineLabel_Internal(tr, label, labelNew));
1995:     PetscCall(DMSetField_Internal(rdm, f, labelNew, obj));
1996:     PetscCall(DMLabelDestroy(&labelNew));
1997:   }
1998:   PetscCall(DMGetNumDS(dm, &Nds));
1999:   for (s = 0; s < Nds; ++s) {
2000:     DMLabel     label, labelNew;
2001:     const char *lname;

2003:     PetscCall(DMGetRegionNumDS(rdm, s, &label, NULL, NULL, NULL));
2004:     if (!label) continue;
2005:     PetscCall(PetscObjectGetName((PetscObject)label, &lname));
2006:     PetscCall(DMLabelCreate(PETSC_COMM_SELF, lname, &labelNew));
2007:     PetscCall(RefineLabel_Internal(tr, label, labelNew));
2008:     PetscCall(DMSetRegionNumDS(rdm, s, labelNew, NULL, NULL, NULL));
2009:     PetscCall(DMLabelDestroy(&labelNew));
2010:   }
2011:   PetscFunctionReturn(PETSC_SUCCESS);
2012: }

2014: static PetscErrorCode DMPlexTransformCreateSF(DMPlexTransform tr, DM rdm)
2015: {
2016:   DM                 dm;
2017:   PetscSF            sf, sfNew;
2018:   PetscInt           numRoots, numLeaves, numLeavesNew = 0, l, m;
2019:   const PetscInt    *localPoints;
2020:   const PetscSFNode *remotePoints;
2021:   PetscInt          *localPointsNew;
2022:   PetscSFNode       *remotePointsNew;
2023:   PetscInt           pStartNew, pEndNew, pNew;
2024:   /* Brute force algorithm */
2025:   PetscSF         rsf;
2026:   PetscSection    s;
2027:   const PetscInt *rootdegree;
2028:   PetscInt       *rootPointsNew, *remoteOffsets;
2029:   PetscInt        numPointsNew, pStart, pEnd, p;

2031:   PetscFunctionBegin;
2032:   PetscCall(DMPlexTransformGetDM(tr, &dm));
2033:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_CreateSF, tr, dm, 0, 0));
2034:   PetscCall(DMPlexGetChart(rdm, &pStartNew, &pEndNew));
2035:   PetscCall(DMGetPointSF(dm, &sf));
2036:   PetscCall(DMGetPointSF(rdm, &sfNew));
2037:   /* Calculate size of new SF */
2038:   PetscCall(PetscSFGetGraph(sf, &numRoots, &numLeaves, &localPoints, &remotePoints));
2039:   if (numRoots < 0) {
2040:     PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_CreateSF, tr, dm, 0, 0));
2041:     PetscFunctionReturn(PETSC_SUCCESS);
2042:   }
2043:   for (l = 0; l < numLeaves; ++l) {
2044:     const PetscInt  p = localPoints[l];
2045:     DMPolytopeType  ct;
2046:     DMPolytopeType *rct;
2047:     PetscInt       *rsize, *rcone, *rornt;
2048:     PetscInt        Nct, n;

2050:     PetscCall(DMPlexGetCellType(dm, p, &ct));
2051:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2052:     for (n = 0; n < Nct; ++n) numLeavesNew += rsize[n];
2053:   }
2054:   /* Send new root point numbers
2055:        It is possible to optimize for regular transforms by sending only the cell type offsets, but it seems a needless complication
2056:   */
2057:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2058:   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &s));
2059:   PetscCall(PetscSectionSetChart(s, pStart, pEnd));
2060:   for (p = pStart; p < pEnd; ++p) {
2061:     DMPolytopeType  ct;
2062:     DMPolytopeType *rct;
2063:     PetscInt       *rsize, *rcone, *rornt;
2064:     PetscInt        Nct, n;

2066:     PetscCall(DMPlexGetCellType(dm, p, &ct));
2067:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2068:     for (n = 0; n < Nct; ++n) PetscCall(PetscSectionAddDof(s, p, rsize[n]));
2069:   }
2070:   PetscCall(PetscSectionSetUp(s));
2071:   PetscCall(PetscSectionGetStorageSize(s, &numPointsNew));
2072:   PetscCall(PetscSFCreateRemoteOffsets(sf, s, s, &remoteOffsets));
2073:   PetscCall(PetscSFCreateSectionSF(sf, s, remoteOffsets, s, &rsf));
2074:   PetscCall(PetscFree(remoteOffsets));
2075:   PetscCall(PetscSFComputeDegreeBegin(sf, &rootdegree));
2076:   PetscCall(PetscSFComputeDegreeEnd(sf, &rootdegree));
2077:   PetscCall(PetscMalloc1(numPointsNew, &rootPointsNew));
2078:   for (p = 0; p < numPointsNew; ++p) rootPointsNew[p] = -1;
2079:   for (p = pStart; p < pEnd; ++p) {
2080:     DMPolytopeType  ct;
2081:     DMPolytopeType *rct;
2082:     PetscInt       *rsize, *rcone, *rornt;
2083:     PetscInt        Nct, n, r, off;

2085:     if (!rootdegree[p - pStart]) continue;
2086:     PetscCall(PetscSectionGetOffset(s, p, &off));
2087:     PetscCall(DMPlexGetCellType(dm, p, &ct));
2088:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2089:     for (n = 0, m = 0; n < Nct; ++n) {
2090:       for (r = 0; r < rsize[n]; ++r, ++m) {
2091:         PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, r, &pNew));
2092:         rootPointsNew[off + m] = pNew;
2093:       }
2094:     }
2095:   }
2096:   PetscCall(PetscSFBcastBegin(rsf, MPIU_INT, rootPointsNew, rootPointsNew, MPI_REPLACE));
2097:   PetscCall(PetscSFBcastEnd(rsf, MPIU_INT, rootPointsNew, rootPointsNew, MPI_REPLACE));
2098:   PetscCall(PetscSFDestroy(&rsf));
2099:   PetscCall(PetscMalloc1(numLeavesNew, &localPointsNew));
2100:   PetscCall(PetscMalloc1(numLeavesNew, &remotePointsNew));
2101:   for (l = 0, m = 0; l < numLeaves; ++l) {
2102:     const PetscInt  p = localPoints[l];
2103:     DMPolytopeType  ct;
2104:     DMPolytopeType *rct;
2105:     PetscInt       *rsize, *rcone, *rornt;
2106:     PetscInt        Nct, n, r, q, off;

2108:     PetscCall(PetscSectionGetOffset(s, p, &off));
2109:     PetscCall(DMPlexGetCellType(dm, p, &ct));
2110:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2111:     for (n = 0, q = 0; n < Nct; ++n) {
2112:       for (r = 0; r < rsize[n]; ++r, ++m, ++q) {
2113:         PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, r, &pNew));
2114:         localPointsNew[m]        = pNew;
2115:         remotePointsNew[m].index = rootPointsNew[off + q];
2116:         remotePointsNew[m].rank  = remotePoints[l].rank;
2117:       }
2118:     }
2119:   }
2120:   PetscCall(PetscSectionDestroy(&s));
2121:   PetscCall(PetscFree(rootPointsNew));
2122:   /* SF needs sorted leaves to correctly calculate Gather */
2123:   {
2124:     PetscSFNode *rp, *rtmp;
2125:     PetscInt    *lp, *idx, *ltmp, i;

2127:     PetscCall(PetscMalloc1(numLeavesNew, &idx));
2128:     PetscCall(PetscMalloc1(numLeavesNew, &lp));
2129:     PetscCall(PetscMalloc1(numLeavesNew, &rp));
2130:     for (i = 0; i < numLeavesNew; ++i) {
2131:       PetscCheck(!(localPointsNew[i] < pStartNew) && !(localPointsNew[i] >= pEndNew), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local SF point %" PetscInt_FMT " (%" PetscInt_FMT ") not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", localPointsNew[i], i, pStartNew, pEndNew);
2132:       idx[i] = i;
2133:     }
2134:     PetscCall(PetscSortIntWithPermutation(numLeavesNew, localPointsNew, idx));
2135:     for (i = 0; i < numLeavesNew; ++i) {
2136:       lp[i] = localPointsNew[idx[i]];
2137:       rp[i] = remotePointsNew[idx[i]];
2138:     }
2139:     ltmp            = localPointsNew;
2140:     localPointsNew  = lp;
2141:     rtmp            = remotePointsNew;
2142:     remotePointsNew = rp;
2143:     PetscCall(PetscFree(idx));
2144:     PetscCall(PetscFree(ltmp));
2145:     PetscCall(PetscFree(rtmp));
2146:   }
2147:   PetscCall(PetscSFSetGraph(sfNew, pEndNew - pStartNew, numLeavesNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
2148:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_CreateSF, tr, dm, 0, 0));
2149:   if (PetscDefined(USE_DEBUG)) {
2150:     PetscInt overlap;

2152:     // Need to set overlap because some transforms put cells in the overlap
2153:     PetscCall(DMPlexGetOverlap(rdm, &overlap));
2154:     PetscCall(DMPlexSetOverlap(rdm, NULL, 1));
2155:     PetscCall(DMPlexCheckPointSF(rdm, sfNew, PETSC_FALSE));
2156:     PetscCall(DMPlexSetOverlap(rdm, NULL, overlap));
2157:   }
2158:   PetscFunctionReturn(PETSC_SUCCESS);
2159: }

2161: /*
2162:   DMPlexCellRefinerMapLocalizedCoordinates - Given a cell of `DMPolytopeType` `ct` with localized coordinates `x`, generate localized coordinates `xr` for subcell `r` of type `rct`.

2164:   Not Collective

2166:   Input Parameters:
2167: + tr  - The `DMPlexTransform`
2168: . ct  - The type of the parent cell
2169: . rct - The type of the produced cell
2170: . r   - The index of the produced cell
2171: - x   - The localized coordinates for the parent cell

2173:   Output Parameter:
2174: . xr  - The localized coordinates for the produced cell

2176:   Level: developer

2178: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPolytopeType`, `DMPlexCellRefinerSetCoordinates()`
2179: */
2180: static PetscErrorCode DMPlexTransformMapLocalizedCoordinates(DMPlexTransform tr, DMPolytopeType ct, DMPolytopeType rct, PetscInt r, const PetscScalar x[], PetscScalar xr[])
2181: {
2182:   PetscFE  fe = NULL;
2183:   PetscInt cdim, v, *subcellV;

2185:   PetscFunctionBegin;
2186:   PetscCall(DMPlexTransformGetCoordinateFE(tr, ct, &fe));
2187:   PetscCall(DMPlexTransformGetSubcellVertices(tr, ct, rct, r, &subcellV));
2188:   PetscCall(PetscFEGetNumComponents(fe, &cdim));
2189:   for (v = 0; v < DMPolytopeTypeGetNumVertices(rct); ++v) PetscCall(PetscFEInterpolate_Static(fe, x, tr->refGeom[ct], subcellV[v], &xr[v * cdim]));
2190:   PetscFunctionReturn(PETSC_SUCCESS);
2191: }

2193: static PetscErrorCode DMPlexTransformSetCoordinates(DMPlexTransform tr, DM rdm)
2194: {
2195:   DM                 dm, cdm, cdmCell, cdmNew, cdmCellNew;
2196:   PetscSection       coordSection, coordSectionNew, coordSectionCell, coordSectionCellNew;
2197:   Vec                coordsLocal, coordsLocalNew, coordsLocalCell = NULL, coordsLocalCellNew;
2198:   const PetscScalar *coords;
2199:   PetscScalar       *coordsNew;
2200:   const PetscReal   *maxCell, *Lstart, *L;
2201:   PetscBool          localized, localizeVertices = PETSC_FALSE, localizeCells = PETSC_FALSE, sparseLocalize;
2202:   PetscInt           dE, dEo, d, cStart, cEnd, c, cStartNew, cEndNew, vStartNew, vEndNew, v, pStart, pEnd, p;

2204:   PetscFunctionBegin;
2205:   // Need to clear the DMField for coordinates
2206:   PetscCall(DMSetCoordinateField(rdm, NULL));
2207:   PetscCall(DMPlexTransformGetDM(tr, &dm));
2208:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_SetCoordinates, tr, dm, 0, 0));
2209:   PetscCall(DMGetCoordinateDM(dm, &cdm));
2210:   PetscCall(DMGetCellCoordinateDM(dm, &cdmCell));
2211:   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
2212:   PetscCall(DMGetSparseLocalize(dm, &sparseLocalize));
2213:   PetscCall(DMSetSparseLocalize(rdm, sparseLocalize));
2214:   PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
2215:   if (localized) {
2216:     /* Localize coordinates of new vertices */
2217:     localizeVertices = PETSC_TRUE;
2218:     /* If we do not have a mechanism for automatically localizing cell coordinates, we need to compute them explicitly for every divided cell */
2219:     if (!maxCell) localizeCells = PETSC_TRUE;
2220:   }
2221:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
2222:   PetscCall(PetscSectionGetFieldComponents(coordSection, 0, &dEo));
2223:   PetscCall(DMGetCoordinateDim(rdm, &dE));
2224:   if (maxCell) {
2225:     PetscReal *LstartNew, *LNew, *maxCellNew;

2227:     PetscCall(PetscMalloc3(dE, &LstartNew, dE, &LNew, dE, &maxCellNew));
2228:     for (d = 0; d < dEo; ++d) {
2229:       LstartNew[d]  = Lstart[d];
2230:       LNew[d]       = L[d];
2231:       maxCellNew[d] = maxCell[d] / tr->redFactor;
2232:     }
2233:     for (d = dEo; d < dE; ++d) {
2234:       LstartNew[d]  = 0.;
2235:       LNew[d]       = -1.;
2236:       maxCellNew[d] = -1.;
2237:     }
2238:     PetscCall(DMSetPeriodicity(rdm, maxCellNew, LstartNew, LNew));
2239:     PetscCall(PetscFree3(LstartNew, LNew, maxCellNew));
2240:   }
2241:   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)rdm), &coordSectionNew));
2242:   PetscCall(PetscSectionSetNumFields(coordSectionNew, 1));
2243:   PetscCall(PetscSectionSetFieldComponents(coordSectionNew, 0, dE));
2244:   PetscCall(DMPlexGetDepthStratum(rdm, 0, &vStartNew, &vEndNew));
2245:   PetscCall(PetscSectionSetChart(coordSectionNew, vStartNew, vEndNew));
2246:   /* Localization should be inherited */
2247:   /*   Stefano calculates parent cells for each new cell for localization */
2248:   /*   Localized cells need coordinates of closure */
2249:   for (v = vStartNew; v < vEndNew; ++v) {
2250:     PetscCall(PetscSectionSetDof(coordSectionNew, v, dE));
2251:     PetscCall(PetscSectionSetFieldDof(coordSectionNew, v, 0, dE));
2252:   }
2253:   PetscCall(PetscSectionSetUp(coordSectionNew));
2254:   PetscCall(DMSetCoordinateSection(rdm, PETSC_DETERMINE, coordSectionNew));

2256:   if (localizeCells) {
2257:     PetscCall(DMGetCoordinateDM(rdm, &cdmNew));
2258:     PetscCall(DMClone(cdmNew, &cdmCellNew));
2259:     PetscCall(DMSetCellCoordinateDM(rdm, cdmCellNew));
2260:     PetscCall(DMDestroy(&cdmCellNew));

2262:     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)rdm), &coordSectionCellNew));
2263:     PetscCall(PetscSectionSetNumFields(coordSectionCellNew, 1));
2264:     PetscCall(PetscSectionSetFieldComponents(coordSectionCellNew, 0, dE));
2265:     PetscCall(DMPlexGetHeightStratum(rdm, 0, &cStartNew, &cEndNew));
2266:     PetscCall(PetscSectionSetChart(coordSectionCellNew, cStartNew, cEndNew));

2268:     PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell));
2269:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2270:     for (c = cStart; c < cEnd; ++c) {
2271:       PetscInt dof;

2273:       PetscCall(PetscSectionGetDof(coordSectionCell, c, &dof));
2274:       if (dof) {
2275:         DMPolytopeType  ct;
2276:         DMPolytopeType *rct;
2277:         PetscInt       *rsize, *rcone, *rornt;
2278:         PetscInt        dim, cNew, Nct, n, r;

2280:         PetscCall(DMPlexGetCellType(dm, c, &ct));
2281:         dim = DMPolytopeTypeGetDim(ct);
2282:         PetscCall(DMPlexTransformCellTransform(tr, ct, c, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2283:         /* This allows for different cell types */
2284:         for (n = 0; n < Nct; ++n) {
2285:           if (dim != DMPolytopeTypeGetDim(rct[n])) continue;
2286:           for (r = 0; r < rsize[n]; ++r) {
2287:             PetscInt *closure = NULL;
2288:             PetscInt  clSize, cl, Nv = 0;

2290:             PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], c, r, &cNew));
2291:             PetscCall(DMPlexGetTransitiveClosure(rdm, cNew, PETSC_TRUE, &clSize, &closure));
2292:             for (cl = 0; cl < clSize * 2; cl += 2) {
2293:               if ((closure[cl] >= vStartNew) && (closure[cl] < vEndNew)) ++Nv;
2294:             }
2295:             PetscCall(DMPlexRestoreTransitiveClosure(rdm, cNew, PETSC_TRUE, &clSize, &closure));
2296:             PetscCall(PetscSectionSetDof(coordSectionCellNew, cNew, Nv * dE));
2297:             PetscCall(PetscSectionSetFieldDof(coordSectionCellNew, cNew, 0, Nv * dE));
2298:           }
2299:         }
2300:       }
2301:     }
2302:     PetscCall(PetscSectionSetUp(coordSectionCellNew));
2303:     PetscCall(DMSetCellCoordinateSection(rdm, PETSC_DETERMINE, coordSectionCellNew));
2304:   }
2305:   PetscCall(DMViewFromOptions(dm, NULL, "-coarse_dm_view"));
2306:   {
2307:     VecType     vtype;
2308:     PetscInt    coordSizeNew, bs;
2309:     const char *name;

2311:     PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
2312:     PetscCall(VecCreate(PETSC_COMM_SELF, &coordsLocalNew));
2313:     PetscCall(PetscSectionGetStorageSize(coordSectionNew, &coordSizeNew));
2314:     PetscCall(VecSetSizes(coordsLocalNew, coordSizeNew, PETSC_DETERMINE));
2315:     PetscCall(PetscObjectGetName((PetscObject)coordsLocal, &name));
2316:     PetscCall(PetscObjectSetName((PetscObject)coordsLocalNew, name));
2317:     PetscCall(VecGetBlockSize(coordsLocal, &bs));
2318:     PetscCall(VecSetBlockSize(coordsLocalNew, dEo == dE ? bs : dE));
2319:     PetscCall(VecGetType(coordsLocal, &vtype));
2320:     PetscCall(VecSetType(coordsLocalNew, vtype));
2321:   }
2322:   PetscCall(VecGetArrayRead(coordsLocal, &coords));
2323:   PetscCall(VecGetArray(coordsLocalNew, &coordsNew));
2324:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2325:   /* First set coordinates for vertices */
2326:   for (p = pStart; p < pEnd; ++p) {
2327:     DMPolytopeType  ct;
2328:     DMPolytopeType *rct;
2329:     PetscInt       *rsize, *rcone, *rornt;
2330:     PetscInt        Nct, n, r;
2331:     PetscBool       hasVertex = PETSC_FALSE;

2333:     PetscCall(DMPlexGetCellType(dm, p, &ct));
2334:     PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2335:     for (n = 0; n < Nct; ++n) {
2336:       if (rct[n] == DM_POLYTOPE_POINT) {
2337:         hasVertex = PETSC_TRUE;
2338:         break;
2339:       }
2340:     }
2341:     if (hasVertex) {
2342:       const PetscScalar *icoords = NULL;
2343:       const PetscScalar *array   = NULL;
2344:       PetscScalar       *pcoords = NULL;
2345:       PetscBool          isDG;
2346:       PetscInt           Nc, Nv, v, d;

2348:       PetscCall(DMPlexGetCellCoordinates(dm, p, &isDG, &Nc, &array, &pcoords));

2350:       icoords = pcoords;
2351:       Nv      = Nc / dEo;
2352:       if (ct != DM_POLYTOPE_POINT) {
2353:         if (localizeVertices && maxCell) {
2354:           PetscScalar anchor[3];

2356:           for (d = 0; d < dEo; ++d) anchor[d] = pcoords[d];
2357:           for (v = 0; v < Nv; ++v) PetscCall(DMLocalizeCoordinate_Internal(dm, dEo, anchor, &pcoords[v * dEo], &pcoords[v * dEo]));
2358:         }
2359:       }
2360:       for (n = 0; n < Nct; ++n) {
2361:         if (rct[n] != DM_POLYTOPE_POINT) continue;
2362:         for (r = 0; r < rsize[n]; ++r) {
2363:           PetscScalar vcoords[3];
2364:           PetscInt    vNew, off;

2366:           PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, r, &vNew));
2367:           PetscCall(PetscSectionGetOffset(coordSectionNew, vNew, &off));
2368:           PetscCall(DMPlexTransformMapCoordinates(tr, ct, rct[n], p, r, Nv, dEo, icoords, vcoords));
2369:           PetscCall(DMSnapToGeomModel(dm, p, dE, vcoords, &coordsNew[off]));
2370:         }
2371:       }
2372:       PetscCall(DMPlexRestoreCellCoordinates(dm, p, &isDG, &Nc, &array, &pcoords));
2373:     }
2374:   }
2375:   PetscCall(VecRestoreArrayRead(coordsLocal, &coords));
2376:   PetscCall(VecRestoreArray(coordsLocalNew, &coordsNew));
2377:   PetscCall(DMSetCoordinatesLocal(rdm, coordsLocalNew));
2378:   PetscCall(VecDestroy(&coordsLocalNew));
2379:   PetscCall(PetscSectionDestroy(&coordSectionNew));
2380:   /* Then set coordinates for cells by localizing */
2381:   if (!localizeCells) PetscCall(DMLocalizeCoordinates(rdm));
2382:   else {
2383:     VecType     vtype;
2384:     PetscInt    coordSizeNew, bs;
2385:     const char *name;

2387:     PetscCall(DMGetCellCoordinatesLocal(dm, &coordsLocalCell));
2388:     PetscCall(VecCreate(PETSC_COMM_SELF, &coordsLocalCellNew));
2389:     PetscCall(PetscSectionGetStorageSize(coordSectionCellNew, &coordSizeNew));
2390:     PetscCall(VecSetSizes(coordsLocalCellNew, coordSizeNew, PETSC_DETERMINE));
2391:     PetscCall(PetscObjectGetName((PetscObject)coordsLocalCell, &name));
2392:     PetscCall(PetscObjectSetName((PetscObject)coordsLocalCellNew, name));
2393:     PetscCall(VecGetBlockSize(coordsLocalCell, &bs));
2394:     PetscCall(VecSetBlockSize(coordsLocalCellNew, dEo == dE ? bs : dE));
2395:     PetscCall(VecGetType(coordsLocalCell, &vtype));
2396:     PetscCall(VecSetType(coordsLocalCellNew, vtype));
2397:     PetscCall(VecGetArrayRead(coordsLocalCell, &coords));
2398:     PetscCall(VecGetArray(coordsLocalCellNew, &coordsNew));

2400:     for (p = pStart; p < pEnd; ++p) {
2401:       DMPolytopeType  ct;
2402:       DMPolytopeType *rct;
2403:       PetscInt       *rsize, *rcone, *rornt;
2404:       PetscInt        dof = 0, Nct, n, r;

2406:       PetscCall(DMPlexGetCellType(dm, p, &ct));
2407:       PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
2408:       if (p >= cStart && p < cEnd) PetscCall(PetscSectionGetDof(coordSectionCell, p, &dof));
2409:       if (dof) {
2410:         const PetscScalar *pcoords;

2412:         PetscCall(DMPlexPointLocalRead(cdmCell, p, coords, &pcoords));
2413:         for (n = 0; n < Nct; ++n) {
2414:           const PetscInt Nr = rsize[n];

2416:           if (DMPolytopeTypeGetDim(ct) != DMPolytopeTypeGetDim(rct[n])) continue;
2417:           for (r = 0; r < Nr; ++r) {
2418:             PetscInt pNew, offNew;

2420:             /* It looks like Stefano and Lisandro are allowing localized coordinates without defining the periodic boundary, which means that
2421:                DMLocalizeCoordinate_Internal() will not work. Localized coordinates will have to have obtained by the affine map of the larger
2422:                cell to the ones it produces. */
2423:             PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, r, &pNew));
2424:             PetscCall(PetscSectionGetOffset(coordSectionCellNew, pNew, &offNew));
2425:             PetscCall(DMPlexTransformMapLocalizedCoordinates(tr, ct, rct[n], r, pcoords, &coordsNew[offNew]));
2426:           }
2427:         }
2428:       }
2429:     }
2430:     PetscCall(VecRestoreArrayRead(coordsLocalCell, &coords));
2431:     PetscCall(VecRestoreArray(coordsLocalCellNew, &coordsNew));
2432:     PetscCall(DMSetCellCoordinatesLocal(rdm, coordsLocalCellNew));
2433:     PetscCall(VecDestroy(&coordsLocalCellNew));
2434:     PetscCall(PetscSectionDestroy(&coordSectionCellNew));
2435:   }
2436:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_SetCoordinates, tr, dm, 0, 0));
2437:   PetscFunctionReturn(PETSC_SUCCESS);
2438: }

2440: /*@
2441:   DMPlexTransformApply - Execute the transformation, producing another `DM`

2443:   Collective

2445:   Input Parameters:
2446: + tr - The `DMPlexTransform` object
2447: - dm - The original `DM`

2449:   Output Parameter:
2450: . trdm - The transformed `DM`

2452:   Level: intermediate

2454:   Options Database Keys:
2455: + -dm_plex_transform_label_match_strata    - Only label points of the same stratum as the producing point
2456: . -dm_plex_transform_label_replica_inc num - Increment for the label value to be multiplied by the replica number
2457: - -dm_plex_transform_active name           - Name for active mesh label

2459: .seealso: [](plex_transform_table), [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransform`, `DMPlexTransformCreate()`, `DMPlexTransformSetDM()`
2460: @*/
2461: PetscErrorCode DMPlexTransformApply(DMPlexTransform tr, DM dm, DM *trdm)
2462: {
2463:   DM                     rdm;
2464:   DMPlexInterpolatedFlag interp;
2465:   PetscInt               pStart, pEnd;

2467:   PetscFunctionBegin;
2470:   PetscAssertPointer(trdm, 3);
2471:   PetscCall(PetscLogEventBegin(DMPLEXTRANSFORM_Apply, tr, dm, 0, 0));
2472:   PetscCall(DMPlexTransformSetDM(tr, dm));

2474:   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &rdm));
2475:   PetscCall(DMSetType(rdm, DMPLEX));
2476:   PetscCall(DMPlexTransformSetDimensions(tr, dm, rdm));
2477:   /* Calculate number of new points of each depth */
2478:   PetscCall(DMPlexIsInterpolatedCollective(dm, &interp));
2479:   PetscCheck(interp == DMPLEX_INTERPOLATED_FULL, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be fully interpolated for regular refinement");
2480:   /* Step 1: Set chart */
2481:   PetscCall(DMPlexTransformGetChart(tr, &pStart, &pEnd));
2482:   PetscCall(DMPlexSetChart(rdm, pStart, pEnd));
2483:   /* Step 2: Set cone/support sizes (automatically stratifies) */
2484:   PetscCall(DMPlexTransformSetConeSizes(tr, rdm));
2485:   /* Step 3: Setup refined DM */
2486:   PetscCall(DMSetUp(rdm));
2487:   /* Step 4: Set cones and supports (automatically symmetrizes) */
2488:   PetscCall(DMPlexTransformSetCones(tr, rdm));
2489:   /* Step 5: Create pointSF */
2490:   PetscCall(DMPlexTransformCreateSF(tr, rdm));
2491:   /* Step 6: Create labels */
2492:   PetscCall(DMPlexTransformCreateLabels(tr, rdm));
2493:   /* Step 7: Set coordinates */
2494:   PetscCall(DMPlexTransformSetCoordinates(tr, rdm));
2495:   //   Do not copy periodicity, which was handled in DMPlexTransformSetCoordinates()
2496:   PetscCall(DMPlexCopy_Internal(dm, PETSC_FALSE, PETSC_TRUE, rdm));
2497:   // If the original DM was configured from options, the transformed DM should be as well
2498:   rdm->setfromoptionscalled = dm->setfromoptionscalled;
2499:   PetscCall(PetscLogEventEnd(DMPLEXTRANSFORM_Apply, tr, dm, 0, 0));
2500:   *trdm = rdm;
2501:   PetscFunctionReturn(PETSC_SUCCESS);
2502: }

2504: PetscErrorCode DMPlexTransformAdaptLabel(DM dm, PETSC_UNUSED Vec metric, DMLabel adaptLabel, PETSC_UNUSED DMLabel rgLabel, DM *rdm)
2505: {
2506:   DMPlexTransform tr;
2507:   DM              cdm, rcdm;
2508:   const char     *prefix;
2509:   PetscBool       save;

2511:   PetscFunctionBegin;
2512:   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
2513:   PetscCall(PetscObjectSetName((PetscObject)tr, "Adapt Label Transform"));
2514:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
2515:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
2516:   PetscCall(DMPlexTransformSetDM(tr, dm));
2517:   PetscCall(DMPlexTransformSetFromOptions(tr));
2518:   if (adaptLabel) PetscCall(DMPlexTransformSetActive(tr, adaptLabel));
2519:   PetscCall(DMPlexTransformSetUp(tr));
2520:   PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
2521:   PetscCall(DMPlexTransformApply(tr, dm, rdm));
2522:   PetscCall(DMCopyDisc(dm, *rdm));
2523:   PetscCall(DMGetCoordinateDM(dm, &cdm));
2524:   PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
2525:   PetscCall(DMCopyDisc(cdm, rcdm));
2526:   PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm));
2527:   PetscCall(DMCopyDisc(dm, *rdm));
2528:   PetscCall(DMPlexGetSaveTransform(dm, &save));
2529:   if (save) PetscCall(DMPlexSetTransform(*rdm, tr));
2530:   PetscCall(DMPlexTransformDestroy(&tr));
2531:   ((DM_Plex *)(*rdm)->data)->useHashLocation = ((DM_Plex *)dm->data)->useHashLocation;
2532:   PetscFunctionReturn(PETSC_SUCCESS);
2533: }