xref: /petsc/src/dm/impls/plex/plexceed.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1f918ec44SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"          I*/
2f918ec44SMatthew G. Knepley 
352e7713aSMatthew G. Knepley static PetscErrorCode DMGetPoints_Private(DM dm, DMLabel domainLabel, PetscInt labelVal, PetscInt height, IS *pointIS)
452e7713aSMatthew G. Knepley {
552e7713aSMatthew G. Knepley   PetscInt depth;
652e7713aSMatthew G. Knepley   DMLabel  depthLabel;
752e7713aSMatthew G. Knepley   IS       depthIS;
852e7713aSMatthew G. Knepley 
952e7713aSMatthew G. Knepley   PetscFunctionBegin;
1052e7713aSMatthew G. Knepley   PetscCall(DMPlexGetDepth(dm, &depth));
1152e7713aSMatthew G. Knepley   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
1252e7713aSMatthew G. Knepley   PetscCall(DMLabelGetStratumIS(depthLabel, depth - height, &depthIS));
1352e7713aSMatthew G. Knepley   if (domainLabel) {
1452e7713aSMatthew G. Knepley     IS domainIS;
1552e7713aSMatthew G. Knepley 
1652e7713aSMatthew G. Knepley     PetscCall(DMLabelGetStratumIS(domainLabel, labelVal, &domainIS));
1752e7713aSMatthew G. Knepley     if (domainIS) { // domainIS is non-empty
1852e7713aSMatthew G. Knepley       PetscCall(ISIntersect(depthIS, domainIS, pointIS));
1952e7713aSMatthew G. Knepley       PetscCall(ISDestroy(&domainIS));
2052e7713aSMatthew G. Knepley     } else { // domainIS is NULL (empty)
2152e7713aSMatthew G. Knepley       *pointIS = NULL;
2252e7713aSMatthew G. Knepley     }
2352e7713aSMatthew G. Knepley     PetscCall(ISDestroy(&depthIS));
2452e7713aSMatthew G. Knepley   } else {
2552e7713aSMatthew G. Knepley     *pointIS = depthIS;
2652e7713aSMatthew G. Knepley   }
2752e7713aSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2852e7713aSMatthew G. Knepley }
2952e7713aSMatthew G. Knepley 
30a2c9b50fSJeremy L Thompson /*@C
3152e7713aSMatthew G. Knepley   DMPlexGetLocalOffsets - Allocate and populate array of local offsets for each cell closure.
3252e7713aSMatthew G. Knepley 
3352e7713aSMatthew G. Knepley   Not collective
34a2c9b50fSJeremy L Thompson 
35a2c9b50fSJeremy L Thompson   Input Parameters:
36a1cb98faSBarry Smith +  dm - The `DMPLEX` object
37a1cb98faSBarry Smith .  domain_label - label for `DMPLEX` domain, or NULL for whole domain
38a1cb98faSBarry Smith .  label_value - Stratum value
39a1cb98faSBarry Smith .  height - Height of target cells in `DMPLEX` topology
40a1cb98faSBarry Smith -  dm_field - Index of `DMPLEX` field
41a2c9b50fSJeremy L Thompson 
42a2c9b50fSJeremy L Thompson   Output Parameters:
43a1cb98faSBarry Smith +  num_cells - Number of local cells
44a1cb98faSBarry Smith .  cell_size - Size of each cell, given by cell_size * num_comp = num_dof
45a1cb98faSBarry Smith .  num_comp - Number of components per dof
46a1cb98faSBarry Smith .  l_size - Size of local vector
47a1cb98faSBarry Smith -  offsets - Allocated offsets array for cells
48a2c9b50fSJeremy L Thompson 
49a2c9b50fSJeremy L Thompson   Level: developer
50a2c9b50fSJeremy L Thompson 
51a1cb98faSBarry Smith   Notes:
52a1cb98faSBarry Smith   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].
53a1cb98faSBarry Smith 
54a1cb98faSBarry Smith    Caller is responsible for freeing the offsets array using `PetscFree()`.
55a1cb98faSBarry Smith 
56*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPlexGetLocalOffsetsSupport()`, `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()`
57a2c9b50fSJeremy L Thompson @*/
58d71ae5a4SJacob Faibussowitsch 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)
59d71ae5a4SJacob Faibussowitsch {
60a2c9b50fSJeremy L Thompson   PetscDS         ds = NULL;
61a2c9b50fSJeremy L Thompson   PetscFE         fe;
62a2c9b50fSJeremy L Thompson   PetscSection    section;
637a17eebaSJed Brown   PetscInt        dim, ds_field = -1;
64a2c9b50fSJeremy L Thompson   PetscInt       *restr_indices;
65a2c9b50fSJeremy L Thompson   const PetscInt *iter_indices;
66a2c9b50fSJeremy L Thompson   IS              iter_is;
67a2c9b50fSJeremy L Thompson 
685f80ce2aSJacob Faibussowitsch   PetscFunctionBeginUser;
695f80ce2aSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
709566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
727a17eebaSJed Brown   {
737a17eebaSJed Brown     IS              field_is;
747a17eebaSJed Brown     const PetscInt *fields;
757a17eebaSJed Brown     PetscInt        num_fields;
765f80ce2aSJacob Faibussowitsch 
7707218a29SMatthew G. Knepley     PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
78a2c9b50fSJeremy L Thompson     // Translate dm_field to ds_field
799566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(field_is, &fields));
809566063dSJacob Faibussowitsch     PetscCall(ISGetSize(field_is, &num_fields));
817a17eebaSJed Brown     for (PetscInt i = 0; i < num_fields; i++) {
827a17eebaSJed Brown       if (dm_field == fields[i]) {
837a17eebaSJed Brown         ds_field = i;
84a2c9b50fSJeremy L Thompson         break;
85a2c9b50fSJeremy L Thompson       }
86a2c9b50fSJeremy L Thompson     }
879566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(field_is, &fields));
88a2c9b50fSJeremy L Thompson   }
8963a3b9bcSJacob Faibussowitsch   PetscCheck(ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
90a2c9b50fSJeremy L Thompson 
9152e7713aSMatthew G. Knepley   PetscCall(DMGetPoints_Private(dm, domain_label, label_value, height, &iter_is));
92a2c9b50fSJeremy L Thompson   if (iter_is) {
939566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(iter_is, num_cells));
949566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(iter_is, &iter_indices));
95a2c9b50fSJeremy L Thompson   } else {
96a2c9b50fSJeremy L Thompson     *num_cells   = 0;
97a2c9b50fSJeremy L Thompson     iter_indices = NULL;
98a2c9b50fSJeremy L Thompson   }
99a2c9b50fSJeremy L Thompson 
100a2c9b50fSJeremy L Thompson   {
101a2c9b50fSJeremy L Thompson     PetscDualSpace dual_space;
102a2c9b50fSJeremy L Thompson     PetscInt       num_dual_basis_vectors;
1035f80ce2aSJacob Faibussowitsch 
1049566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
1059566063dSJacob Faibussowitsch     PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
1067c48043bSMatthew G. Knepley     PetscCheck(fe, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " is invalid for DG coordinates", height);
1079566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
1089566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
1099566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetNumComponents(dual_space, num_comp));
11063a3b9bcSJacob Faibussowitsch     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);
111a2c9b50fSJeremy L Thompson     *cell_size = num_dual_basis_vectors / *num_comp;
112a2c9b50fSJeremy L Thompson   }
113a2c9b50fSJeremy L Thompson   PetscInt restr_size = (*num_cells) * (*cell_size);
1149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(restr_size, &restr_indices));
115a2c9b50fSJeremy L Thompson   PetscInt cell_offset = 0;
116a2c9b50fSJeremy L Thompson 
11752e7713aSMatthew G. Knepley   PetscInt P = (PetscInt)PetscPowReal(*cell_size, 1.0 / (dim - height));
118a2c9b50fSJeremy L Thompson   for (PetscInt p = 0; p < *num_cells; p++) {
119a2c9b50fSJeremy L Thompson     PetscBool flip = PETSC_FALSE;
120a2c9b50fSJeremy L Thompson     PetscInt  c    = iter_indices[p];
121a2c9b50fSJeremy L Thompson     PetscInt  num_indices, *indices;
122a2c9b50fSJeremy L Thompson     PetscInt  field_offsets[17]; // max number of fields plus 1
1239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
124a2c9b50fSJeremy L Thompson     if (height > 0) {
125a2c9b50fSJeremy L Thompson       PetscInt        num_cells_support, num_faces, start = -1;
126a2c9b50fSJeremy L Thompson       const PetscInt *orients, *faces, *cells;
1279566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, c, &cells));
1289566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, c, &num_cells_support));
1295f80ce2aSJacob Faibussowitsch       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);
1309566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, cells[0], &faces));
1319566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, cells[0], &num_faces));
1329371c9d4SSatish Balay       for (PetscInt i = 0; i < num_faces; i++) {
1339371c9d4SSatish Balay         if (faces[i] == c) start = i;
1349371c9d4SSatish Balay       }
1355f80ce2aSJacob Faibussowitsch       PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c);
1369566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeOrientation(dm, cells[0], &orients));
137a2c9b50fSJeremy L Thompson       if (orients[start] < 0) flip = PETSC_TRUE;
138a2c9b50fSJeremy L Thompson     }
139a2c9b50fSJeremy L Thompson 
140a2c9b50fSJeremy L Thompson     for (PetscInt i = 0; i < *cell_size; i++) {
141a2c9b50fSJeremy L Thompson       PetscInt ii = i;
142a2c9b50fSJeremy L Thompson       if (flip) {
143a2c9b50fSJeremy L Thompson         if (*cell_size == P) ii = *cell_size - 1 - i;
144a2c9b50fSJeremy L Thompson         else if (*cell_size == P * P) {
145a2c9b50fSJeremy L Thompson           PetscInt row = i / P, col = i % P;
146a2c9b50fSJeremy L Thompson           ii = row + col * P;
14763a3b9bcSJacob Faibussowitsch         } 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);
148a2c9b50fSJeremy L Thompson       }
149a2c9b50fSJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
150a2c9b50fSJeremy L Thompson       PetscInt loc                 = indices[field_offsets[dm_field] + ii * (*num_comp)];
151a2c9b50fSJeremy L Thompson       restr_indices[cell_offset++] = loc >= 0 ? loc : -(loc + 1);
152a2c9b50fSJeremy L Thompson     }
1539566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
154a2c9b50fSJeremy L Thompson   }
15563a3b9bcSJacob Faibussowitsch   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);
1569566063dSJacob Faibussowitsch   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
1579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iter_is));
158a2c9b50fSJeremy L Thompson 
159a2c9b50fSJeremy L Thompson   *offsets = restr_indices;
1609566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, l_size));
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162a2c9b50fSJeremy L Thompson }
163a2c9b50fSJeremy L Thompson 
16452e7713aSMatthew G. Knepley /*@C
16552e7713aSMatthew G. Knepley   DMPlexGetLocalOffsetsSupport - Allocate and populate arrays of local offsets for each face support.
16652e7713aSMatthew G. Knepley 
16752e7713aSMatthew G. Knepley   Not collective
16852e7713aSMatthew G. Knepley 
16952e7713aSMatthew G. Knepley   Input Parameters:
17052e7713aSMatthew G. Knepley +  dm - The `DMPLEX` object
17152e7713aSMatthew G. Knepley .  domain_label - label for `DMPLEX` domain, or NULL for whole domain
17252e7713aSMatthew G. Knepley -  label_value - Stratum value
17352e7713aSMatthew G. Knepley 
17452e7713aSMatthew G. Knepley   Output Parameters:
17552e7713aSMatthew G. Knepley +  num_faces - Number of local, non-boundary faces
17652e7713aSMatthew G. Knepley .  num_comp - Number of components per dof
17752e7713aSMatthew G. Knepley .  l_size - Size of local vector
17852e7713aSMatthew G. Knepley .  offsetsNeg - Allocated offsets array for cells on the inward normal side of each face
17952e7713aSMatthew G. Knepley -  offsetsPos - Allocated offsets array for cells on the outward normal side of each face
18052e7713aSMatthew G. Knepley 
18152e7713aSMatthew G. Knepley   Level: developer
18252e7713aSMatthew G. Knepley 
18352e7713aSMatthew G. Knepley   Notes:
18452e7713aSMatthew G. Knepley   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].
18552e7713aSMatthew G. Knepley 
18652e7713aSMatthew G. Knepley    Caller is responsible for freeing the offsets array using `PetscFree()`.
18752e7713aSMatthew G. Knepley 
188*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPlexGetLocalOffsets()`, `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()`
18952e7713aSMatthew G. Knepley @*/
19052e7713aSMatthew G. Knepley PetscErrorCode DMPlexGetLocalOffsetsSupport(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt *num_faces, PetscInt *num_comp, PetscInt *l_size, PetscInt **offsetsNeg, PetscInt **offsetsPos)
19152e7713aSMatthew G. Knepley {
19252e7713aSMatthew G. Knepley   PetscDS         ds = NULL;
19352e7713aSMatthew G. Knepley   PetscFV         fv;
19452e7713aSMatthew G. Knepley   PetscSection    section;
19552e7713aSMatthew G. Knepley   PetscInt        dim, height = 1, dm_field = 0, ds_field = 0, Nf, NfInt = 0, cell_size, restr_size;
19652e7713aSMatthew G. Knepley   PetscInt       *restr_indices_neg, *restr_indices_pos;
19752e7713aSMatthew G. Knepley   const PetscInt *iter_indices;
19852e7713aSMatthew G. Knepley   IS              iter_is;
19952e7713aSMatthew G. Knepley 
20052e7713aSMatthew G. Knepley   PetscFunctionBeginUser;
20152e7713aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
20252e7713aSMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
20352e7713aSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
20407218a29SMatthew G. Knepley   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
20552e7713aSMatthew G. Knepley 
20652e7713aSMatthew G. Knepley   PetscCall(DMGetPoints_Private(dm, domain_label, label_value, height, &iter_is));
20752e7713aSMatthew G. Knepley   if (iter_is) {
20852e7713aSMatthew G. Knepley     PetscCall(ISGetIndices(iter_is, &iter_indices));
20952e7713aSMatthew G. Knepley     PetscCall(ISGetLocalSize(iter_is, &Nf));
21052e7713aSMatthew G. Knepley     for (PetscInt p = 0, Ns; p < Nf; ++p) {
21152e7713aSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dm, iter_indices[p], &Ns));
21252e7713aSMatthew G. Knepley       if (Ns == 2) ++NfInt;
21352e7713aSMatthew G. Knepley     }
21452e7713aSMatthew G. Knepley     *num_faces = NfInt;
21552e7713aSMatthew G. Knepley   } else {
21652e7713aSMatthew G. Knepley     *num_faces   = 0;
21752e7713aSMatthew G. Knepley     iter_indices = NULL;
21852e7713aSMatthew G. Knepley   }
21952e7713aSMatthew G. Knepley 
22052e7713aSMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fv));
22152e7713aSMatthew G. Knepley   PetscCall(PetscFVGetNumComponents(fv, num_comp));
22252e7713aSMatthew G. Knepley   cell_size  = *num_comp;
22352e7713aSMatthew G. Knepley   restr_size = (*num_faces) * cell_size;
22452e7713aSMatthew G. Knepley   PetscCall(PetscMalloc1(restr_size, &restr_indices_neg));
22552e7713aSMatthew G. Knepley   PetscCall(PetscMalloc1(restr_size, &restr_indices_pos));
22652e7713aSMatthew G. Knepley   PetscInt face_offset_neg = 0, face_offset_pos = 0;
22752e7713aSMatthew G. Knepley 
22852e7713aSMatthew G. Knepley   for (PetscInt p = 0; p < Nf; ++p) {
22952e7713aSMatthew G. Knepley     const PetscInt  face = iter_indices[p];
23052e7713aSMatthew G. Knepley     PetscInt        num_indices, *indices;
23152e7713aSMatthew G. Knepley     PetscInt        field_offsets[17]; // max number of fields plus 1
23252e7713aSMatthew G. Knepley     const PetscInt *supp;
23352e7713aSMatthew G. Knepley     PetscInt        Ns;
23452e7713aSMatthew G. Knepley 
23552e7713aSMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, face, &supp));
23652e7713aSMatthew G. Knepley     PetscCall(DMPlexGetSupportSize(dm, face, &Ns));
23752e7713aSMatthew G. Knepley     // Ignore boundary faces
23852e7713aSMatthew G. Knepley     //   TODO check for face on parallel boundary
23952e7713aSMatthew G. Knepley     if (Ns == 2) {
24052e7713aSMatthew G. Knepley       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
24152e7713aSMatthew G. Knepley       PetscCall(DMPlexGetClosureIndices(dm, section, section, supp[0], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
24252e7713aSMatthew G. Knepley       for (PetscInt i = 0; i < cell_size; i++) {
24352e7713aSMatthew G. Knepley         const PetscInt loc                   = indices[field_offsets[dm_field] + i * (*num_comp)];
24452e7713aSMatthew G. Knepley         restr_indices_neg[face_offset_neg++] = loc >= 0 ? loc : -(loc + 1);
24552e7713aSMatthew G. Knepley       }
24652e7713aSMatthew G. Knepley       PetscCall(DMPlexRestoreClosureIndices(dm, section, section, supp[0], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
24752e7713aSMatthew G. Knepley       PetscCall(DMPlexGetClosureIndices(dm, section, section, supp[1], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
24852e7713aSMatthew G. Knepley       for (PetscInt i = 0; i < cell_size; i++) {
24952e7713aSMatthew G. Knepley         const PetscInt loc                   = indices[field_offsets[dm_field] + i * (*num_comp)];
25052e7713aSMatthew G. Knepley         restr_indices_pos[face_offset_pos++] = loc >= 0 ? loc : -(loc + 1);
25152e7713aSMatthew G. Knepley       }
25252e7713aSMatthew G. Knepley       PetscCall(DMPlexRestoreClosureIndices(dm, section, section, supp[1], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
25352e7713aSMatthew G. Knepley     }
25452e7713aSMatthew G. Knepley   }
25552e7713aSMatthew G. Knepley   PetscCheck(face_offset_neg == restr_size, PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, neg offsets array of shape (%" PetscInt_FMT ", %" PetscInt_FMT ") initialized for %" PetscInt_FMT " nodes", *num_faces, cell_size, face_offset_neg);
25652e7713aSMatthew G. Knepley   PetscCheck(face_offset_pos == restr_size, PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, pos offsets array of shape (%" PetscInt_FMT ", %" PetscInt_FMT ") initialized for %" PetscInt_FMT " nodes", *num_faces, cell_size, face_offset_pos);
25752e7713aSMatthew G. Knepley   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
25852e7713aSMatthew G. Knepley   PetscCall(ISDestroy(&iter_is));
25952e7713aSMatthew G. Knepley 
26052e7713aSMatthew G. Knepley   *offsetsNeg = restr_indices_neg;
26152e7713aSMatthew G. Knepley   *offsetsPos = restr_indices_pos;
26252e7713aSMatthew G. Knepley   PetscCall(PetscSectionGetStorageSize(section, l_size));
26352e7713aSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
26452e7713aSMatthew G. Knepley }
26552e7713aSMatthew G. Knepley 
266a2c9b50fSJeremy L Thompson #if defined(PETSC_HAVE_LIBCEED)
267f918ec44SMatthew G. Knepley   #include <petscdmplexceed.h>
268f918ec44SMatthew G. Knepley 
269a2c9b50fSJeremy L Thompson /*@C
270a2c9b50fSJeremy L Thompson   DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector)
271a2c9b50fSJeremy L Thompson 
272a2c9b50fSJeremy L Thompson   Input Parameters:
273a1cb98faSBarry Smith +  dm - The `DMPLEX` object
274a1cb98faSBarry Smith .  domain_label - label for `DMPLEX` domain, or NULL for the whole domain
275a1cb98faSBarry Smith .  label_value - Stratum value
276a1cb98faSBarry Smith .  height - Height of target cells in `DMPLEX` topology
277a1cb98faSBarry Smith -  dm_field - Index of `DMPLEX` field
278a2c9b50fSJeremy L Thompson 
279a1cb98faSBarry Smith   Output Parameter:
280a1cb98faSBarry Smith .  ERestrict - libCEED restriction from local vector to to the cells
281a2c9b50fSJeremy L Thompson 
282a2c9b50fSJeremy L Thompson   Level: developer
283a1cb98faSBarry Smith 
284*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetLocalOffsets()`
285a2c9b50fSJeremy L Thompson @*/
286d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict)
287d71ae5a4SJacob Faibussowitsch {
288f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
289f918ec44SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
29021066fe8SJed Brown   PetscValidPointer(ERestrict, 6);
291f918ec44SMatthew G. Knepley   if (!dm->ceedERestrict) {
292a2c9b50fSJeremy L Thompson     PetscInt            num_cells, cell_size, num_comp, lvec_size, *restr_indices;
293a2c9b50fSJeremy L Thompson     CeedElemRestriction elem_restr;
294f918ec44SMatthew G. Knepley     Ceed                ceed;
295f918ec44SMatthew G. Knepley 
2969566063dSJacob Faibussowitsch     PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices));
2979566063dSJacob Faibussowitsch     PetscCall(DMGetCeed(dm, &ceed));
2989566063dSJacob Faibussowitsch     PetscCallCEED(CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, restr_indices, &elem_restr));
2999566063dSJacob Faibussowitsch     PetscCall(PetscFree(restr_indices));
300a2c9b50fSJeremy L Thompson     dm->ceedERestrict = elem_restr;
301f918ec44SMatthew G. Knepley   }
302f918ec44SMatthew G. Knepley   *ERestrict = dm->ceedERestrict;
3033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
304f918ec44SMatthew G. Knepley }
305f918ec44SMatthew G. Knepley 
306f918ec44SMatthew G. Knepley #endif
307