Actual source code: plexceed.c

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

  3: PetscErrorCode DMGetPoints_Internal(DM dm, DMLabel domainLabel, PetscInt labelVal, PetscInt height, IS *pointIS)
  4: {
  5:   PetscInt depth;
  6:   DMLabel  depthLabel;
  7:   IS       depthIS;

  9:   PetscFunctionBegin;
 10:   PetscCall(DMPlexGetDepth(dm, &depth));
 11:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
 12:   PetscCall(DMLabelGetStratumIS(depthLabel, depth - height, &depthIS));
 13:   if (domainLabel) {
 14:     IS domainIS;

 16:     PetscCall(DMLabelGetStratumIS(domainLabel, labelVal, &domainIS));
 17:     if (domainIS) { // domainIS is non-empty
 18:       PetscCall(ISIntersect(depthIS, domainIS, pointIS));
 19:       PetscCall(ISDestroy(&domainIS));
 20:     } else { // domainIS is NULL (empty)
 21:       *pointIS = NULL;
 22:     }
 23:     PetscCall(ISDestroy(&depthIS));
 24:   } else {
 25:     *pointIS = depthIS;
 26:   }
 27:   PetscFunctionReturn(PETSC_SUCCESS);
 28: }

 30: /*@C
 31:   DMPlexGetLocalOffsets - Allocate and populate array of local offsets for each cell closure.

 33:   Not collective

 35:   Input Parameters:
 36: + dm           - The `DMPLEX` object
 37: . domain_label - label for `DMPLEX` domain, or NULL for whole domain
 38: . label_value  - Stratum value
 39: . height       - Height of target cells in `DMPLEX` topology
 40: - dm_field     - Index of `DMPLEX` field

 42:   Output Parameters:
 43: + num_cells - Number of local cells
 44: . cell_size - Size of each cell, given by cell_size * num_comp = num_dof
 45: . num_comp  - Number of components per dof
 46: . l_size    - Size of local vector
 47: - offsets   - Allocated offsets array for cells

 49:   Level: developer

 51:   Notes:
 52:   Allocate and populate array of shape [num_cells, cell_size] defining offsets for each value (cell, node) for local vector of the `DMPLEX` field. All offsets are in the range [0, l_size - 1].

 54:   Caller is responsible for freeing the offsets array using `PetscFree()`.

 56: .seealso: [](ch_unstructured), `DMPlexGetLocalOffsetsSupport()`, `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()`
 57: @*/
 58: PetscErrorCode DMPlexGetLocalOffsets(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, PetscInt *num_cells, PetscInt *cell_size, PetscInt *num_comp, PetscInt *l_size, PetscInt **offsets)
 59: {
 60:   PetscDS         ds = NULL;
 61:   PetscFE         fe;
 62:   PetscSection    section;
 63:   PetscInt        dim, ds_field = -1;
 64:   PetscInt       *restr_indices;
 65:   const PetscInt *iter_indices;
 66:   IS              iter_is;

 68:   PetscFunctionBeginUser;
 70:   PetscCall(PetscLogEventBegin(DMPLEX_GetLocalOffsets, dm, 0, 0, 0));
 71:   PetscCall(DMGetDimension(dm, &dim));
 72:   PetscCall(DMGetLocalSection(dm, &section));
 73:   PetscCall(PetscSectionGetStorageSize(section, l_size));
 74:   {
 75:     IS              field_is;
 76:     const PetscInt *fields;
 77:     PetscInt        num_fields;

 79:     PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
 80:     // Translate dm_field to ds_field
 81:     PetscCall(ISGetIndices(field_is, &fields));
 82:     PetscCall(ISGetSize(field_is, &num_fields));
 83:     for (PetscInt i = 0; i < num_fields; i++) {
 84:       if (dm_field == fields[i]) {
 85:         ds_field = i;
 86:         break;
 87:       }
 88:     }
 89:     PetscCall(ISRestoreIndices(field_is, &fields));
 90:   }
 91:   PetscCheck(ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);

 93:   PetscCall(DMGetPoints_Internal(dm, domain_label, label_value, height, &iter_is));
 94:   if (iter_is) {
 95:     PetscCall(ISGetLocalSize(iter_is, num_cells));
 96:     PetscCall(ISGetIndices(iter_is, &iter_indices));
 97:   } else {
 98:     *num_cells   = 0;
 99:     iter_indices = NULL;
100:   }

102:   {
103:     PetscDualSpace dual_space;
104:     PetscInt       num_dual_basis_vectors;

106:     PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
107:     PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
108:     PetscCheck(fe, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " is invalid for DG discretizations", height);
109:     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
110:     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
111:     PetscCall(PetscDualSpaceGetNumComponents(dual_space, num_comp));
112:     PetscCheck(num_dual_basis_vectors % *num_comp == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for number of dual basis vectors %" PetscInt_FMT " not divisible by %" PetscInt_FMT " components", num_dual_basis_vectors, *num_comp);
113:     *cell_size = num_dual_basis_vectors / *num_comp;
114:   }
115:   PetscInt restr_size = (*num_cells) * (*cell_size);
116:   PetscCall(PetscMalloc1(restr_size, &restr_indices));
117:   PetscInt cell_offset = 0;

119:   PetscInt P = dim - height ? (PetscInt)PetscPowReal(*cell_size, 1.0 / (dim - height)) : 0;
120:   for (PetscInt p = 0; p < *num_cells; p++) {
121:     PetscBool flip = PETSC_FALSE;
122:     PetscInt  c    = iter_indices[p];
123:     PetscInt  num_indices, *indices;
124:     PetscInt  field_offsets[17]; // max number of fields plus 1
125:     PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
126:     if (height > 0) {
127:       PetscInt        num_cells_support, num_faces, start = -1;
128:       const PetscInt *orients, *faces, *cells;
129:       PetscCall(DMPlexGetSupport(dm, c, &cells));
130:       PetscCall(DMPlexGetSupportSize(dm, c, &num_cells_support));
131:       PetscCheck(num_cells_support == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Expected one cell in support of exterior face, but got %" PetscInt_FMT " cells", num_cells_support);
132:       PetscCall(DMPlexGetCone(dm, cells[0], &faces));
133:       PetscCall(DMPlexGetConeSize(dm, cells[0], &num_faces));
134:       for (PetscInt i = 0; i < num_faces; i++) {
135:         if (faces[i] == c) start = i;
136:       }
137:       PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c);
138:       PetscCall(DMPlexGetConeOrientation(dm, cells[0], &orients));
139:       if (orients[start] < 0) flip = PETSC_TRUE;
140:     }

142:     for (PetscInt i = 0; i < *cell_size; i++) {
143:       PetscInt ii = i;
144:       if (flip) {
145:         if (*cell_size == P) ii = *cell_size - 1 - i;
146:         else if (*cell_size == P * P) {
147:           PetscInt row = i / P, col = i % P;
148:           ii = row + col * P;
149:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for flipping point with cell size %" PetscInt_FMT " != P (%" PetscInt_FMT ") or P^2", *cell_size, P);
150:       }
151:       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
152:       PetscInt loc                 = indices[field_offsets[dm_field] + ii * (*num_comp)];
153:       loc                          = loc < 0 ? -(loc + 1) : loc;
154:       restr_indices[cell_offset++] = loc;
155:       PetscCheck(loc >= 0 && loc < *l_size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Location %" PetscInt_FMT " not in [0, %" PetscInt_FMT ") local vector", loc, *l_size);
156:     }
157:     PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
158:   }
159:   PetscCheck(cell_offset == restr_size, PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, offsets array of shape (%" PetscInt_FMT ", %" PetscInt_FMT ") initialized for %" PetscInt_FMT " nodes", *num_cells, *cell_size, cell_offset);
160:   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
161:   PetscCall(ISDestroy(&iter_is));

163:   *offsets = restr_indices;
164:   PetscCall(PetscLogEventEnd(DMPLEX_GetLocalOffsets, dm, 0, 0, 0));
165:   PetscFunctionReturn(PETSC_SUCCESS);
166: }

168: /*@C
169:   DMPlexGetLocalOffsetsSupport - Allocate and populate arrays of local offsets for each face support.

171:   Not collective

173:   Input Parameters:
174: + dm           - The `DMPLEX` object
175: . domain_label - label for `DMPLEX` domain, or NULL for whole domain
176: - label_value  - Stratum value

178:   Output Parameters:
179: + num_faces  - Number of local, non-boundary faces
180: . num_comp   - Number of components per dof
181: . l_size     - Size of local vector
182: . offsetsNeg - Allocated offsets array for cells on the inward normal side of each face
183: - offsetsPos - Allocated offsets array for cells on the outward normal side of each face

185:   Level: developer

187:   Notes:
188:   Allocate and populate array of shape [num_cells, num_comp] defining offsets for each cell for local vector of the `DMPLEX` field. All offsets are in the range [0, l_size - 1].

190:   Caller is responsible for freeing the offsets array using `PetscFree()`.

192: .seealso: [](ch_unstructured), `DMPlexGetLocalOffsets()`, `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()`
193: @*/
194: PetscErrorCode DMPlexGetLocalOffsetsSupport(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt *num_faces, PetscInt *num_comp, PetscInt *l_size, PetscInt **offsetsNeg, PetscInt **offsetsPos)
195: {
196:   PetscDS         ds = NULL;
197:   PetscFV         fv;
198:   PetscSection    section;
199:   PetscInt        dim, height = 1, dm_field = 0, ds_field = 0, Nf, NfInt = 0, Nc;
200:   PetscInt       *restr_indices_neg, *restr_indices_pos;
201:   const PetscInt *iter_indices;
202:   IS              iter_is;

204:   PetscFunctionBeginUser;
206:   PetscCall(DMGetDimension(dm, &dim));
207:   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
208:   PetscCall(DMGetLocalSection(dm, &section));
209:   PetscCall(PetscSectionGetStorageSize(section, l_size));

211:   PetscCall(DMGetPoints_Internal(dm, domain_label, label_value, height, &iter_is));
212:   if (iter_is) {
213:     PetscCall(ISGetIndices(iter_is, &iter_indices));
214:     PetscCall(ISGetLocalSize(iter_is, &Nf));
215:     for (PetscInt p = 0, Ns; p < Nf; ++p) {
216:       PetscCall(DMPlexGetSupportSize(dm, iter_indices[p], &Ns));
217:       if (Ns == 2) ++NfInt;
218:     }
219:     *num_faces = NfInt;
220:   } else {
221:     *num_faces   = 0;
222:     iter_indices = NULL;
223:   }

225:   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fv));
226:   PetscCall(PetscFVGetNumComponents(fv, &Nc));
227:   PetscCall(PetscMalloc1(NfInt, &restr_indices_neg));
228:   PetscCall(PetscMalloc1(NfInt, &restr_indices_pos));
229:   PetscInt face_offset_neg = 0, face_offset_pos = 0;

231:   for (PetscInt p = 0; p < Nf; ++p) {
232:     const PetscInt  face = iter_indices[p];
233:     PetscInt        num_indices, *indices;
234:     PetscInt        field_offsets[17]; // max number of fields plus 1
235:     const PetscInt *supp;
236:     PetscInt        Ns, loc;

238:     PetscCall(DMPlexGetSupport(dm, face, &supp));
239:     PetscCall(DMPlexGetSupportSize(dm, face, &Ns));
240:     // Ignore boundary faces
241:     //   TODO check for face on parallel boundary
242:     if (Ns == 2) {
243:       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
244:       PetscCall(DMPlexGetClosureIndices(dm, section, section, supp[0], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
245:       PetscCheck(num_indices == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of closure indices %" PetscInt_FMT " != %" PetscInt_FMT " number of FV components", num_indices, Nc);
246:       loc                                  = indices[field_offsets[dm_field]];
247:       loc                                  = loc < 0 ? -(loc + 1) : loc;
248:       restr_indices_neg[face_offset_neg++] = loc;
249:       PetscCheck(loc >= 0 && loc + Nc <= *l_size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Location %" PetscInt_FMT " + Nc not in [0, %" PetscInt_FMT ") local vector", loc, *l_size);
250:       PetscCall(DMPlexRestoreClosureIndices(dm, section, section, supp[0], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
251:       PetscCall(DMPlexGetClosureIndices(dm, section, section, supp[1], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
252:       PetscCheck(num_indices == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of closure indices %" PetscInt_FMT " != %" PetscInt_FMT " number of FV components", num_indices, Nc);
253:       loc                                  = indices[field_offsets[dm_field]];
254:       loc                                  = loc < 0 ? -(loc + 1) : loc;
255:       restr_indices_pos[face_offset_pos++] = loc;
256:       PetscCheck(loc >= 0 && loc + Nc <= *l_size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Location %" PetscInt_FMT " + Nc not in [0, %" PetscInt_FMT ") local vector", loc, *l_size);
257:       PetscCall(DMPlexRestoreClosureIndices(dm, section, section, supp[1], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
258:     }
259:   }
260:   PetscCheck(face_offset_neg == NfInt, PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, neg offsets array of shape (%" PetscInt_FMT ") initialized for %" PetscInt_FMT " nodes", NfInt, face_offset_neg);
261:   PetscCheck(face_offset_pos == NfInt, PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, pos offsets array of shape (%" PetscInt_FMT ") initialized for %" PetscInt_FMT " nodes", NfInt, face_offset_pos);
262:   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
263:   PetscCall(ISDestroy(&iter_is));

265:   *num_comp   = Nc;
266:   *offsetsNeg = restr_indices_neg;
267:   *offsetsPos = restr_indices_pos;
268:   PetscFunctionReturn(PETSC_SUCCESS);
269: }

271: #if defined(PETSC_HAVE_LIBCEED)
272: #include <petscdmplexceed.h>

274: // Consumes the input petsc_indices and provides the output ceed_indices; no-copy when the sizes match.
275: static PetscErrorCode PetscIntArrayIntoCeedInt_Private(PetscInt length, PetscInt max_bound, const char *max_bound_name, PetscInt **petsc_indices, CeedInt **ceed_indices)
276: {
277:   PetscFunctionBegin;
278:   if (length) PetscAssertPointer(petsc_indices, 3);
279:   PetscAssertPointer(ceed_indices, 4);
280:   #if defined(PETSC_USE_64BIT_INDICES)
281:   PetscCheck(max_bound <= PETSC_INT32_MAX, PETSC_COMM_SELF, PETSC_ERR_SUP, "%s %" PetscInt_FMT " does not fit in int32_t", max_bound_name, max_bound);
282:   {
283:     CeedInt *ceed_ind;
284:     PetscCall(PetscMalloc1(length, &ceed_ind));
285:     for (PetscInt i = 0; i < length; i++) ceed_ind[i] = (*petsc_indices)[i];
286:     *ceed_indices = ceed_ind;
287:     PetscCall(PetscFree(*petsc_indices));
288:   }
289:   #else
290:   *ceed_indices  = *petsc_indices;
291:   *petsc_indices = NULL;
292:   #endif
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

296: /*@C
297:   DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector)

299:   Input Parameters:
300: +  dm - The `DMPLEX` object
301: .  domain_label - label for `DMPLEX` domain, or NULL for the whole domain
302: .  label_value - Stratum value
303: .  height - Height of target cells in `DMPLEX` topology
304: -  dm_field - Index of `DMPLEX` field

306:   Output Parameter:
307: .  ERestrict - libCEED restriction from local vector to the cells

309:   Level: developer

311: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetLocalOffsets()`
312: @*/
313: PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict)
314: {
315:   PetscFunctionBeginUser;
317:   PetscAssertPointer(ERestrict, 6);
318:   if (!dm->ceedERestrict) {
319:     PetscInt            num_cells, cell_size, num_comp, lvec_size, *restr_indices;
320:     CeedElemRestriction elem_restr;
321:     Ceed                ceed;

323:     PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices));
324:     PetscCall(DMGetCeed(dm, &ceed));
325:     {
326:       CeedInt *ceed_indices;
327:       PetscCall(PetscIntArrayIntoCeedInt_Private(num_cells * cell_size, lvec_size, "lvec_size", &restr_indices, &ceed_indices));
328:       PetscCallCEED(CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, ceed_indices, &elem_restr));
329:       PetscCall(PetscFree(ceed_indices));
330:     }
331:     dm->ceedERestrict = elem_restr;
332:   }
333:   *ERestrict = dm->ceedERestrict;
334:   PetscFunctionReturn(PETSC_SUCCESS);
335: }

337: PetscErrorCode DMPlexCreateCeedRestrictionFVM(DM dm, CeedElemRestriction *erL, CeedElemRestriction *erR)
338: {
339:   Ceed      ceed;
340:   PetscInt *offL, *offR;
341:   PetscInt  num_faces, num_comp, lvec_size;

343:   PetscFunctionBeginUser;
345:   PetscAssertPointer(erL, 2);
346:   PetscAssertPointer(erR, 3);
347:   PetscCall(DMGetCeed(dm, &ceed));
348:   PetscCall(DMPlexGetLocalOffsetsSupport(dm, NULL, 0, &num_faces, &num_comp, &lvec_size, &offL, &offR));
349:   {
350:     CeedInt *ceed_off;
351:     PetscCall(PetscIntArrayIntoCeedInt_Private(num_faces * 1, lvec_size, "lvec_size", &offL, &ceed_off));
352:     PetscCallCEED(CeedElemRestrictionCreate(ceed, num_faces, 1, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, ceed_off, erL));
353:     PetscCall(PetscFree(ceed_off));
354:     PetscCall(PetscIntArrayIntoCeedInt_Private(num_faces * 1, lvec_size, "lvec_size", &offR, &ceed_off));
355:     PetscCallCEED(CeedElemRestrictionCreate(ceed, num_faces, 1, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, ceed_off, erR));
356:     PetscCall(PetscFree(ceed_off));
357:   }
358:   PetscFunctionReturn(PETSC_SUCCESS);
359: }

361: // TODO DMPlexComputeGeometryFVM() also computes centroids and minimum radius
362: // TODO DMPlexComputeGeometryFVM() flips normal to match support orientation
363: // This function computes area-weights normals
364: PetscErrorCode DMPlexCeedComputeGeometryFVM(DM dm, CeedVector qd)
365: {
366:   DMLabel         domain_label = NULL;
367:   PetscInt        label_value = 0, height = 1, Nf, NfInt = 0, cdim;
368:   const PetscInt *iter_indices;
369:   IS              iter_is;
370:   CeedScalar     *qdata;

372:   PetscFunctionBegin;
373:   PetscCall(DMGetCoordinateDim(dm, &cdim));
374:   PetscCall(DMGetPoints_Internal(dm, domain_label, label_value, height, &iter_is));
375:   if (iter_is) {
376:     PetscCall(ISGetIndices(iter_is, &iter_indices));
377:     PetscCall(ISGetLocalSize(iter_is, &Nf));
378:     for (PetscInt p = 0, Ns; p < Nf; ++p) {
379:       PetscCall(DMPlexGetSupportSize(dm, iter_indices[p], &Ns));
380:       if (Ns == 2) ++NfInt;
381:     }
382:   } else {
383:     iter_indices = NULL;
384:   }

386:   PetscCallCEED(CeedVectorSetValue(qd, 0.));
387:   PetscCallCEED(CeedVectorGetArray(qd, CEED_MEM_HOST, &qdata));
388:   for (PetscInt p = 0, off = 0; p < Nf; ++p) {
389:     const PetscInt  face = iter_indices[p];
390:     const PetscInt *supp;
391:     PetscInt        suppSize;

393:     PetscCall(DMPlexGetSupport(dm, face, &supp));
394:     PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
395:     // Ignore boundary faces
396:     //   TODO check for face on parallel boundary
397:     if (suppSize == 2) {
398:       DMPolytopeType ct;
399:       PetscReal      area, fcentroid[3], centroids[2][3];

401:       PetscCall(DMPlexComputeCellGeometryFVM(dm, face, &area, fcentroid, &qdata[off]));
402:       for (PetscInt d = 0; d < cdim; ++d) qdata[off + d] *= area;
403:       off += cdim;
404:       for (PetscInt s = 0; s < suppSize; ++s) {
405:         PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
406:         if (ct == DM_POLYTOPE_FV_GHOST) continue;
407:         PetscCall(DMPlexComputeCellGeometryFVM(dm, supp[s], &qdata[off + s], centroids[s], NULL));
408:       }
409:       // Give FV ghosts the same volume as the opposite cell
410:       for (PetscInt s = 0; s < suppSize; ++s) {
411:         PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
412:         if (ct != DM_POLYTOPE_FV_GHOST) continue;
413:         qdata[off + s] = qdata[off + (1 - s)];
414:         for (PetscInt d = 0; d < cdim; ++d) centroids[s][d] = fcentroid[d];
415:       }
416:       // Flip normal orientation if necessary to match ordering in support
417:       {
418:         CeedScalar *normal = &qdata[off - cdim];
419:         PetscReal   l[3], r[3], v[3];

421:         PetscCall(DMLocalizeCoordinateReal_Internal(dm, cdim, fcentroid, centroids[0], l));
422:         PetscCall(DMLocalizeCoordinateReal_Internal(dm, cdim, fcentroid, centroids[1], r));
423:         DMPlex_WaxpyD_Internal(cdim, -1, l, r, v);
424:         if (DMPlex_DotRealD_Internal(cdim, normal, v) < 0) {
425:           for (PetscInt d = 0; d < cdim; ++d) normal[d] = -normal[d];
426:         }
427:         if (DMPlex_DotRealD_Internal(cdim, normal, v) <= 0) {
428:           PetscCheck(cdim != 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " could not be fixed, normal (%g,%g) v (%g,%g)", face, (double)normal[0], (double)normal[1], (double)v[0], (double)v[1]);
429:           PetscCheck(cdim != 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", face, (double)normal[0], (double)normal[1], (double)normal[2], (double)v[0], (double)v[1], (double)v[2]);
430:           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " could not be fixed", face);
431:         }
432:       }
433:       off += suppSize;
434:     }
435:   }
436:   PetscCallCEED(CeedVectorRestoreArray(qd, &qdata));
437:   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
438:   PetscCall(ISDestroy(&iter_is));
439:   PetscFunctionReturn(PETSC_SUCCESS);
440: }

442: #endif