xref: /petsc/src/dm/impls/plex/plexceed.c (revision 354b609cae3b349426bee2f9f20157b4134dc6bc)
1f918ec44SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"          I*/
2f918ec44SMatthew G. Knepley 
3*354b609cSMatthew G. Knepley PetscErrorCode DMGetPoints_Internal(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 
561cc06b55SBarry 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);
70708be2fdSJed Brown   PetscCall(PetscLogEventBegin(DMPLEX_GetLocalOffsets, dm, 0, 0, 0));
719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
725962854dSMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
735962854dSMatthew G. Knepley   PetscCall(PetscSectionGetStorageSize(section, l_size));
747a17eebaSJed Brown   {
757a17eebaSJed Brown     IS              field_is;
767a17eebaSJed Brown     const PetscInt *fields;
777a17eebaSJed Brown     PetscInt        num_fields;
785f80ce2aSJacob Faibussowitsch 
7907218a29SMatthew G. Knepley     PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
80a2c9b50fSJeremy L Thompson     // Translate dm_field to ds_field
819566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(field_is, &fields));
829566063dSJacob Faibussowitsch     PetscCall(ISGetSize(field_is, &num_fields));
837a17eebaSJed Brown     for (PetscInt i = 0; i < num_fields; i++) {
847a17eebaSJed Brown       if (dm_field == fields[i]) {
857a17eebaSJed Brown         ds_field = i;
86a2c9b50fSJeremy L Thompson         break;
87a2c9b50fSJeremy L Thompson       }
88a2c9b50fSJeremy L Thompson     }
899566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(field_is, &fields));
90a2c9b50fSJeremy L Thompson   }
9163a3b9bcSJacob Faibussowitsch   PetscCheck(ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
92a2c9b50fSJeremy L Thompson 
93*354b609cSMatthew G. Knepley   PetscCall(DMGetPoints_Internal(dm, domain_label, label_value, height, &iter_is));
94a2c9b50fSJeremy L Thompson   if (iter_is) {
959566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(iter_is, num_cells));
969566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(iter_is, &iter_indices));
97a2c9b50fSJeremy L Thompson   } else {
98a2c9b50fSJeremy L Thompson     *num_cells   = 0;
99a2c9b50fSJeremy L Thompson     iter_indices = NULL;
100a2c9b50fSJeremy L Thompson   }
101a2c9b50fSJeremy L Thompson 
102a2c9b50fSJeremy L Thompson   {
103a2c9b50fSJeremy L Thompson     PetscDualSpace dual_space;
104a2c9b50fSJeremy L Thompson     PetscInt       num_dual_basis_vectors;
1055f80ce2aSJacob Faibussowitsch 
1069566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
1079566063dSJacob Faibussowitsch     PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
1085f82726aSMatthew G. Knepley     PetscCheck(fe, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " is invalid for DG discretizations", height);
1099566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
1109566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
1119566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetNumComponents(dual_space, num_comp));
11263a3b9bcSJacob 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);
113a2c9b50fSJeremy L Thompson     *cell_size = num_dual_basis_vectors / *num_comp;
114a2c9b50fSJeremy L Thompson   }
115a2c9b50fSJeremy L Thompson   PetscInt restr_size = (*num_cells) * (*cell_size);
1169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(restr_size, &restr_indices));
117a2c9b50fSJeremy L Thompson   PetscInt cell_offset = 0;
118a2c9b50fSJeremy L Thompson 
1195f82726aSMatthew G. Knepley   PetscInt P = dim - height ? (PetscInt)PetscPowReal(*cell_size, 1.0 / (dim - height)) : 0;
120a2c9b50fSJeremy L Thompson   for (PetscInt p = 0; p < *num_cells; p++) {
121a2c9b50fSJeremy L Thompson     PetscBool flip = PETSC_FALSE;
122a2c9b50fSJeremy L Thompson     PetscInt  c    = iter_indices[p];
123a2c9b50fSJeremy L Thompson     PetscInt  num_indices, *indices;
124a2c9b50fSJeremy L Thompson     PetscInt  field_offsets[17]; // max number of fields plus 1
1259566063dSJacob Faibussowitsch     PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
126a2c9b50fSJeremy L Thompson     if (height > 0) {
127a2c9b50fSJeremy L Thompson       PetscInt        num_cells_support, num_faces, start = -1;
128a2c9b50fSJeremy L Thompson       const PetscInt *orients, *faces, *cells;
1299566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, c, &cells));
1309566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, c, &num_cells_support));
1315f80ce2aSJacob 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);
1329566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, cells[0], &faces));
1339566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, cells[0], &num_faces));
1349371c9d4SSatish Balay       for (PetscInt i = 0; i < num_faces; i++) {
1359371c9d4SSatish Balay         if (faces[i] == c) start = i;
1369371c9d4SSatish Balay       }
1375f80ce2aSJacob Faibussowitsch       PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c);
1389566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeOrientation(dm, cells[0], &orients));
139a2c9b50fSJeremy L Thompson       if (orients[start] < 0) flip = PETSC_TRUE;
140a2c9b50fSJeremy L Thompson     }
141a2c9b50fSJeremy L Thompson 
142a2c9b50fSJeremy L Thompson     for (PetscInt i = 0; i < *cell_size; i++) {
143a2c9b50fSJeremy L Thompson       PetscInt ii = i;
144a2c9b50fSJeremy L Thompson       if (flip) {
145a2c9b50fSJeremy L Thompson         if (*cell_size == P) ii = *cell_size - 1 - i;
146a2c9b50fSJeremy L Thompson         else if (*cell_size == P * P) {
147a2c9b50fSJeremy L Thompson           PetscInt row = i / P, col = i % P;
148a2c9b50fSJeremy L Thompson           ii = row + col * P;
14963a3b9bcSJacob 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);
150a2c9b50fSJeremy L Thompson       }
151a2c9b50fSJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
152a2c9b50fSJeremy L Thompson       PetscInt loc                 = indices[field_offsets[dm_field] + ii * (*num_comp)];
1535962854dSMatthew G. Knepley       loc                          = loc < 0 ? -(loc + 1) : loc;
1545962854dSMatthew G. Knepley       restr_indices[cell_offset++] = loc;
1555962854dSMatthew G. Knepley       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);
156a2c9b50fSJeremy L Thompson     }
1579566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
158a2c9b50fSJeremy L Thompson   }
15963a3b9bcSJacob 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);
1609566063dSJacob Faibussowitsch   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
1619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iter_is));
162a2c9b50fSJeremy L Thompson 
163a2c9b50fSJeremy L Thompson   *offsets = restr_indices;
164708be2fdSJed Brown   PetscCall(PetscLogEventEnd(DMPLEX_GetLocalOffsets, dm, 0, 0, 0));
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166a2c9b50fSJeremy L Thompson }
167a2c9b50fSJeremy L Thompson 
16852e7713aSMatthew G. Knepley /*@C
16952e7713aSMatthew G. Knepley   DMPlexGetLocalOffsetsSupport - Allocate and populate arrays of local offsets for each face support.
17052e7713aSMatthew G. Knepley 
17152e7713aSMatthew G. Knepley   Not collective
17252e7713aSMatthew G. Knepley 
17352e7713aSMatthew G. Knepley   Input Parameters:
17452e7713aSMatthew G. Knepley + dm           - The `DMPLEX` object
17552e7713aSMatthew G. Knepley . domain_label - label for `DMPLEX` domain, or NULL for whole domain
17652e7713aSMatthew G. Knepley - label_value  - Stratum value
17752e7713aSMatthew G. Knepley 
17852e7713aSMatthew G. Knepley   Output Parameters:
17952e7713aSMatthew G. Knepley + num_faces  - Number of local, non-boundary faces
18052e7713aSMatthew G. Knepley . num_comp   - Number of components per dof
18152e7713aSMatthew G. Knepley . l_size     - Size of local vector
18252e7713aSMatthew G. Knepley . offsetsNeg - Allocated offsets array for cells on the inward normal side of each face
18352e7713aSMatthew G. Knepley - offsetsPos - Allocated offsets array for cells on the outward normal side of each face
18452e7713aSMatthew G. Knepley 
18552e7713aSMatthew G. Knepley   Level: developer
18652e7713aSMatthew G. Knepley 
18752e7713aSMatthew G. Knepley   Notes:
18852e7713aSMatthew 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].
18952e7713aSMatthew G. Knepley 
19052e7713aSMatthew G. Knepley   Caller is responsible for freeing the offsets array using `PetscFree()`.
19152e7713aSMatthew G. Knepley 
1921cc06b55SBarry Smith .seealso: [](ch_unstructured), `DMPlexGetLocalOffsets()`, `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()`
19352e7713aSMatthew G. Knepley @*/
19452e7713aSMatthew 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)
19552e7713aSMatthew G. Knepley {
19652e7713aSMatthew G. Knepley   PetscDS         ds = NULL;
19752e7713aSMatthew G. Knepley   PetscFV         fv;
19852e7713aSMatthew G. Knepley   PetscSection    section;
1995962854dSMatthew G. Knepley   PetscInt        dim, height = 1, dm_field = 0, ds_field = 0, Nf, NfInt = 0, Nc;
20052e7713aSMatthew G. Knepley   PetscInt       *restr_indices_neg, *restr_indices_pos;
20152e7713aSMatthew G. Knepley   const PetscInt *iter_indices;
20252e7713aSMatthew G. Knepley   IS              iter_is;
20352e7713aSMatthew G. Knepley 
20452e7713aSMatthew G. Knepley   PetscFunctionBeginUser;
20552e7713aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
20652e7713aSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
20707218a29SMatthew G. Knepley   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
2085962854dSMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
2095962854dSMatthew G. Knepley   PetscCall(PetscSectionGetStorageSize(section, l_size));
21052e7713aSMatthew G. Knepley 
211*354b609cSMatthew G. Knepley   PetscCall(DMGetPoints_Internal(dm, domain_label, label_value, height, &iter_is));
21252e7713aSMatthew G. Knepley   if (iter_is) {
21352e7713aSMatthew G. Knepley     PetscCall(ISGetIndices(iter_is, &iter_indices));
21452e7713aSMatthew G. Knepley     PetscCall(ISGetLocalSize(iter_is, &Nf));
21552e7713aSMatthew G. Knepley     for (PetscInt p = 0, Ns; p < Nf; ++p) {
21652e7713aSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dm, iter_indices[p], &Ns));
21752e7713aSMatthew G. Knepley       if (Ns == 2) ++NfInt;
21852e7713aSMatthew G. Knepley     }
21952e7713aSMatthew G. Knepley     *num_faces = NfInt;
22052e7713aSMatthew G. Knepley   } else {
22152e7713aSMatthew G. Knepley     *num_faces   = 0;
22252e7713aSMatthew G. Knepley     iter_indices = NULL;
22352e7713aSMatthew G. Knepley   }
22452e7713aSMatthew G. Knepley 
22552e7713aSMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fv));
2265962854dSMatthew G. Knepley   PetscCall(PetscFVGetNumComponents(fv, &Nc));
22762d9cbd2SMatthew G. Knepley   PetscCall(PetscMalloc1(NfInt, &restr_indices_neg));
22862d9cbd2SMatthew G. Knepley   PetscCall(PetscMalloc1(NfInt, &restr_indices_pos));
22952e7713aSMatthew G. Knepley   PetscInt face_offset_neg = 0, face_offset_pos = 0;
23052e7713aSMatthew G. Knepley 
23152e7713aSMatthew G. Knepley   for (PetscInt p = 0; p < Nf; ++p) {
23252e7713aSMatthew G. Knepley     const PetscInt  face = iter_indices[p];
23352e7713aSMatthew G. Knepley     PetscInt        num_indices, *indices;
23452e7713aSMatthew G. Knepley     PetscInt        field_offsets[17]; // max number of fields plus 1
23552e7713aSMatthew G. Knepley     const PetscInt *supp;
2365962854dSMatthew G. Knepley     PetscInt        Ns, loc;
23752e7713aSMatthew G. Knepley 
23852e7713aSMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, face, &supp));
23952e7713aSMatthew G. Knepley     PetscCall(DMPlexGetSupportSize(dm, face, &Ns));
24052e7713aSMatthew G. Knepley     // Ignore boundary faces
24152e7713aSMatthew G. Knepley     //   TODO check for face on parallel boundary
24252e7713aSMatthew G. Knepley     if (Ns == 2) {
24352e7713aSMatthew G. Knepley       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
24452e7713aSMatthew G. Knepley       PetscCall(DMPlexGetClosureIndices(dm, section, section, supp[0], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
2455962854dSMatthew G. Knepley       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);
24662d9cbd2SMatthew G. Knepley       loc                                  = indices[field_offsets[dm_field]];
2475962854dSMatthew G. Knepley       loc                                  = loc < 0 ? -(loc + 1) : loc;
2485962854dSMatthew G. Knepley       restr_indices_neg[face_offset_neg++] = loc;
24962d9cbd2SMatthew G. Knepley       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);
25052e7713aSMatthew G. Knepley       PetscCall(DMPlexRestoreClosureIndices(dm, section, section, supp[0], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
25152e7713aSMatthew G. Knepley       PetscCall(DMPlexGetClosureIndices(dm, section, section, supp[1], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
2525962854dSMatthew G. Knepley       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);
25362d9cbd2SMatthew G. Knepley       loc                                  = indices[field_offsets[dm_field]];
2545962854dSMatthew G. Knepley       loc                                  = loc < 0 ? -(loc + 1) : loc;
2555962854dSMatthew G. Knepley       restr_indices_pos[face_offset_pos++] = loc;
25662d9cbd2SMatthew G. Knepley       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);
25752e7713aSMatthew G. Knepley       PetscCall(DMPlexRestoreClosureIndices(dm, section, section, supp[1], PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
25852e7713aSMatthew G. Knepley     }
25952e7713aSMatthew G. Knepley   }
26062d9cbd2SMatthew G. Knepley   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);
26162d9cbd2SMatthew G. Knepley   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);
26252e7713aSMatthew G. Knepley   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
26352e7713aSMatthew G. Knepley   PetscCall(ISDestroy(&iter_is));
26452e7713aSMatthew G. Knepley 
2655962854dSMatthew G. Knepley   *num_comp   = Nc;
26652e7713aSMatthew G. Knepley   *offsetsNeg = restr_indices_neg;
26752e7713aSMatthew G. Knepley   *offsetsPos = restr_indices_pos;
26852e7713aSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
26952e7713aSMatthew G. Knepley }
27052e7713aSMatthew G. Knepley 
271a2c9b50fSJeremy L Thompson #if defined(PETSC_HAVE_LIBCEED)
272f918ec44SMatthew G. Knepley   #include <petscdmplexceed.h>
273f918ec44SMatthew G. Knepley 
27485f1dce8SJed Brown // Consumes the input petsc_indices and provides the output ceed_indices; no-copy when the sizes match.
27585f1dce8SJed Brown static PetscErrorCode PetscIntArrayIntoCeedInt_Private(PetscInt length, PetscInt max_bound, const char *max_bound_name, PetscInt **petsc_indices, CeedInt **ceed_indices)
27685f1dce8SJed Brown {
27785f1dce8SJed Brown   PetscFunctionBegin;
27885f1dce8SJed Brown   if (length) PetscAssertPointer(petsc_indices, 3);
27985f1dce8SJed Brown   PetscAssertPointer(ceed_indices, 4);
28085f1dce8SJed Brown   #if defined(PETSC_USE_64BIT_INDICES)
28185f1dce8SJed Brown   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);
28285f1dce8SJed Brown   {
28385f1dce8SJed Brown     CeedInt *ceed_ind;
28485f1dce8SJed Brown     PetscCall(PetscMalloc1(length, &ceed_ind));
28585f1dce8SJed Brown     for (PetscInt i = 0; i < length; i++) ceed_ind[i] = (*petsc_indices)[i];
28685f1dce8SJed Brown     *ceed_indices = ceed_ind;
28785f1dce8SJed Brown     PetscCall(PetscFree(*petsc_indices));
28885f1dce8SJed Brown   }
28985f1dce8SJed Brown   #else
29085f1dce8SJed Brown   *ceed_indices  = *petsc_indices;
29185f1dce8SJed Brown   *petsc_indices = NULL;
29285f1dce8SJed Brown   #endif
29385f1dce8SJed Brown   PetscFunctionReturn(PETSC_SUCCESS);
29485f1dce8SJed Brown }
29585f1dce8SJed Brown 
296a2c9b50fSJeremy L Thompson /*@C
297a2c9b50fSJeremy L Thompson   DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector)
298a2c9b50fSJeremy L Thompson 
299a2c9b50fSJeremy L Thompson   Input Parameters:
300a1cb98faSBarry Smith +  dm - The `DMPLEX` object
301a1cb98faSBarry Smith .  domain_label - label for `DMPLEX` domain, or NULL for the whole domain
302a1cb98faSBarry Smith .  label_value - Stratum value
303a1cb98faSBarry Smith .  height - Height of target cells in `DMPLEX` topology
304a1cb98faSBarry Smith -  dm_field - Index of `DMPLEX` field
305a2c9b50fSJeremy L Thompson 
306a1cb98faSBarry Smith   Output Parameter:
307a1cb98faSBarry Smith .  ERestrict - libCEED restriction from local vector to to the cells
308a2c9b50fSJeremy L Thompson 
309a2c9b50fSJeremy L Thompson   Level: developer
310a1cb98faSBarry Smith 
3111cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetLocalOffsets()`
312a2c9b50fSJeremy L Thompson @*/
313d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict)
314d71ae5a4SJacob Faibussowitsch {
315f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
316f918ec44SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3174f572ea9SToby Isaac   PetscAssertPointer(ERestrict, 6);
318f918ec44SMatthew G. Knepley   if (!dm->ceedERestrict) {
319a2c9b50fSJeremy L Thompson     PetscInt            num_cells, cell_size, num_comp, lvec_size, *restr_indices;
320a2c9b50fSJeremy L Thompson     CeedElemRestriction elem_restr;
321f918ec44SMatthew G. Knepley     Ceed                ceed;
322f918ec44SMatthew G. Knepley 
3239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices));
3249566063dSJacob Faibussowitsch     PetscCall(DMGetCeed(dm, &ceed));
32585f1dce8SJed Brown     {
32685f1dce8SJed Brown       CeedInt *ceed_indices;
32785f1dce8SJed Brown       PetscCall(PetscIntArrayIntoCeedInt_Private(num_cells * cell_size, lvec_size, "lvec_size", &restr_indices, &ceed_indices));
32885f1dce8SJed Brown       PetscCallCEED(CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, ceed_indices, &elem_restr));
32985f1dce8SJed Brown       PetscCall(PetscFree(ceed_indices));
33085f1dce8SJed Brown     }
331a2c9b50fSJeremy L Thompson     dm->ceedERestrict = elem_restr;
332f918ec44SMatthew G. Knepley   }
333f918ec44SMatthew G. Knepley   *ERestrict = dm->ceedERestrict;
3343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
335f918ec44SMatthew G. Knepley }
336f918ec44SMatthew G. Knepley 
3375962854dSMatthew G. Knepley PetscErrorCode DMPlexCreateCeedRestrictionFVM(DM dm, CeedElemRestriction *erL, CeedElemRestriction *erR)
3385962854dSMatthew G. Knepley {
3395962854dSMatthew G. Knepley   Ceed      ceed;
3405962854dSMatthew G. Knepley   PetscInt *offL, *offR;
3415962854dSMatthew G. Knepley   PetscInt  num_faces, num_comp, lvec_size;
3425962854dSMatthew G. Knepley 
3435962854dSMatthew G. Knepley   PetscFunctionBeginUser;
3445962854dSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3455962854dSMatthew G. Knepley   PetscAssertPointer(erL, 2);
3465962854dSMatthew G. Knepley   PetscAssertPointer(erR, 3);
3475962854dSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
3485962854dSMatthew G. Knepley   PetscCall(DMPlexGetLocalOffsetsSupport(dm, NULL, 0, &num_faces, &num_comp, &lvec_size, &offL, &offR));
34985f1dce8SJed Brown   {
35085f1dce8SJed Brown     CeedInt *ceed_off;
35185f1dce8SJed Brown     PetscCall(PetscIntArrayIntoCeedInt_Private(num_faces * 1, lvec_size, "lvec_size", &offL, &ceed_off));
35285f1dce8SJed Brown     PetscCallCEED(CeedElemRestrictionCreate(ceed, num_faces, 1, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, ceed_off, erL));
35385f1dce8SJed Brown     PetscCall(PetscFree(ceed_off));
35485f1dce8SJed Brown     PetscCall(PetscIntArrayIntoCeedInt_Private(num_faces * 1, lvec_size, "lvec_size", &offR, &ceed_off));
35585f1dce8SJed Brown     PetscCallCEED(CeedElemRestrictionCreate(ceed, num_faces, 1, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, ceed_off, erR));
35685f1dce8SJed Brown     PetscCall(PetscFree(ceed_off));
35785f1dce8SJed Brown   }
3585962854dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3595962854dSMatthew G. Knepley }
3605962854dSMatthew G. Knepley 
3615962854dSMatthew G. Knepley // TODO DMPlexComputeGeometryFVM() also computes centroids and minimum radius
3625962854dSMatthew G. Knepley // TODO DMPlexComputeGeometryFVM() flips normal to match support orientation
3635962854dSMatthew G. Knepley // This function computes area-weights normals
3645962854dSMatthew G. Knepley PetscErrorCode DMPlexCeedComputeGeometryFVM(DM dm, CeedVector qd)
3655962854dSMatthew G. Knepley {
3665962854dSMatthew G. Knepley   DMLabel         domain_label = NULL;
3675962854dSMatthew G. Knepley   PetscInt        label_value = 0, height = 1, Nf, NfInt = 0, cdim;
3685962854dSMatthew G. Knepley   const PetscInt *iter_indices;
3695962854dSMatthew G. Knepley   IS              iter_is;
3705962854dSMatthew G. Knepley   CeedScalar     *qdata;
3715962854dSMatthew G. Knepley 
3725962854dSMatthew G. Knepley   PetscFunctionBegin;
3735962854dSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
374*354b609cSMatthew G. Knepley   PetscCall(DMGetPoints_Internal(dm, domain_label, label_value, height, &iter_is));
3755962854dSMatthew G. Knepley   if (iter_is) {
3765962854dSMatthew G. Knepley     PetscCall(ISGetIndices(iter_is, &iter_indices));
3775962854dSMatthew G. Knepley     PetscCall(ISGetLocalSize(iter_is, &Nf));
3785962854dSMatthew G. Knepley     for (PetscInt p = 0, Ns; p < Nf; ++p) {
3795962854dSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dm, iter_indices[p], &Ns));
3805962854dSMatthew G. Knepley       if (Ns == 2) ++NfInt;
3815962854dSMatthew G. Knepley     }
3825962854dSMatthew G. Knepley   } else {
3835962854dSMatthew G. Knepley     iter_indices = NULL;
3845962854dSMatthew G. Knepley   }
3855962854dSMatthew G. Knepley 
3865962854dSMatthew G. Knepley   PetscCallCEED(CeedVectorSetValue(qd, 0.));
3875962854dSMatthew G. Knepley   PetscCallCEED(CeedVectorGetArray(qd, CEED_MEM_HOST, &qdata));
3885962854dSMatthew G. Knepley   for (PetscInt p = 0, off = 0; p < Nf; ++p) {
3895962854dSMatthew G. Knepley     const PetscInt  face = iter_indices[p];
3905962854dSMatthew G. Knepley     const PetscInt *supp;
3915962854dSMatthew G. Knepley     PetscInt        suppSize;
3925962854dSMatthew G. Knepley 
3935962854dSMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, face, &supp));
3945962854dSMatthew G. Knepley     PetscCall(DMPlexGetSupportSize(dm, face, &suppSize));
3955962854dSMatthew G. Knepley     // Ignore boundary faces
3965962854dSMatthew G. Knepley     //   TODO check for face on parallel boundary
3975962854dSMatthew G. Knepley     if (suppSize == 2) {
3985962854dSMatthew G. Knepley       DMPolytopeType ct;
3995962854dSMatthew G. Knepley       PetscReal      area, fcentroid[3], centroids[2][3];
4005962854dSMatthew G. Knepley 
4015962854dSMatthew G. Knepley       PetscCall(DMPlexComputeCellGeometryFVM(dm, face, &area, fcentroid, &qdata[off]));
4025962854dSMatthew G. Knepley       for (PetscInt d = 0; d < cdim; ++d) qdata[off + d] *= area;
4035962854dSMatthew G. Knepley       off += cdim;
4045962854dSMatthew G. Knepley       for (PetscInt s = 0; s < suppSize; ++s) {
4055962854dSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
4065962854dSMatthew G. Knepley         if (ct == DM_POLYTOPE_FV_GHOST) continue;
4075962854dSMatthew G. Knepley         PetscCall(DMPlexComputeCellGeometryFVM(dm, supp[s], &qdata[off + s], centroids[s], NULL));
4085962854dSMatthew G. Knepley       }
4095962854dSMatthew G. Knepley       // Give FV ghosts the same volume as the opposite cell
4105962854dSMatthew G. Knepley       for (PetscInt s = 0; s < suppSize; ++s) {
4115962854dSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
4125962854dSMatthew G. Knepley         if (ct != DM_POLYTOPE_FV_GHOST) continue;
4135962854dSMatthew G. Knepley         qdata[off + s] = qdata[off + (1 - s)];
4145962854dSMatthew G. Knepley         for (PetscInt d = 0; d < cdim; ++d) centroids[s][d] = fcentroid[d];
4155962854dSMatthew G. Knepley       }
4165962854dSMatthew G. Knepley       // Flip normal orientation if necessary to match ordering in support
4175962854dSMatthew G. Knepley       {
4185962854dSMatthew G. Knepley         CeedScalar *normal = &qdata[off - cdim];
4195962854dSMatthew G. Knepley         PetscReal   l[3], r[3], v[3];
4205962854dSMatthew G. Knepley 
4215962854dSMatthew G. Knepley         PetscCall(DMLocalizeCoordinateReal_Internal(dm, cdim, fcentroid, centroids[0], l));
4225962854dSMatthew G. Knepley         PetscCall(DMLocalizeCoordinateReal_Internal(dm, cdim, fcentroid, centroids[1], r));
4235962854dSMatthew G. Knepley         DMPlex_WaxpyD_Internal(cdim, -1, l, r, v);
4245962854dSMatthew G. Knepley         if (DMPlex_DotRealD_Internal(cdim, normal, v) < 0) {
4255962854dSMatthew G. Knepley           for (PetscInt d = 0; d < cdim; ++d) normal[d] = -normal[d];
4265962854dSMatthew G. Knepley         }
4275962854dSMatthew G. Knepley         if (DMPlex_DotRealD_Internal(cdim, normal, v) <= 0) {
4285962854dSMatthew G. Knepley           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]);
4295962854dSMatthew G. Knepley           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]);
4305962854dSMatthew G. Knepley           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " could not be fixed", face);
4315962854dSMatthew G. Knepley         }
4325962854dSMatthew G. Knepley       }
4335962854dSMatthew G. Knepley       off += suppSize;
4345962854dSMatthew G. Knepley     }
4355962854dSMatthew G. Knepley   }
4365962854dSMatthew G. Knepley   PetscCallCEED(CeedVectorRestoreArray(qd, &qdata));
4375962854dSMatthew G. Knepley   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
4385962854dSMatthew G. Knepley   PetscCall(ISDestroy(&iter_is));
4395962854dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
4405962854dSMatthew G. Knepley }
4415962854dSMatthew G. Knepley 
442f918ec44SMatthew G. Knepley #endif
443