xref: /petsc/src/dm/impls/plex/plexceed.c (revision 7c48043b109ee9f116d1c17e77d6369f7a4a2b2e)
1f918ec44SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>           /*I      "petscdmplex.h"          I*/
2f918ec44SMatthew G. Knepley 
3a2c9b50fSJeremy L Thompson /*@C
4a2c9b50fSJeremy L Thompson   DMPlexGetLocalOffsets - Allocate and populate array of local offsets.
5a2c9b50fSJeremy L Thompson 
6a2c9b50fSJeremy L Thompson   Input Parameters:
7a2c9b50fSJeremy L Thompson   dm - The DMPlex object
8a2c9b50fSJeremy L Thompson   domain_label - label for DMPlex domain, or NULL for whole domain
9a2c9b50fSJeremy L Thompson   label_value - Stratum value
10a2c9b50fSJeremy L Thompson   height - Height of target cells in DMPlex topology
11a2c9b50fSJeremy L Thompson   dm_field - Index of DMPlex field
12a2c9b50fSJeremy L Thompson 
13a2c9b50fSJeremy L Thompson   Output Parameters:
14a2c9b50fSJeremy L Thompson   num_cells - Number of local cells
15a2c9b50fSJeremy L Thompson   cell_size - Size of each cell, given by cell_size * num_comp = num_dof
16a2c9b50fSJeremy L Thompson   num_comp - Number of components per dof
17a2c9b50fSJeremy L Thompson   l_size - Size of local vector
18a2c9b50fSJeremy L Thompson   offsets - Allocated offsets array for cells
19a2c9b50fSJeremy L Thompson 
20a2c9b50fSJeremy L Thompson   Notes: 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]. Caller is responsible for freeing the offsets array using PetscFree().
21a2c9b50fSJeremy L Thompson 
22a2c9b50fSJeremy L Thompson   Level: developer
23a2c9b50fSJeremy L Thompson 
24db781477SPatrick Sanan .seealso: `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()`
25a2c9b50fSJeremy L Thompson @*/
26a2c9b50fSJeremy L Thompson 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)
27a2c9b50fSJeremy L Thompson {
28a2c9b50fSJeremy L Thompson   PetscDS         ds = NULL;
29a2c9b50fSJeremy L Thompson   PetscFE         fe;
30a2c9b50fSJeremy L Thompson   PetscSection    section;
317a17eebaSJed Brown   PetscInt        dim, ds_field = -1;
32a2c9b50fSJeremy L Thompson   PetscInt       *restr_indices;
33a2c9b50fSJeremy L Thompson   const PetscInt *iter_indices;
34a2c9b50fSJeremy L Thompson   IS              iter_is;
35a2c9b50fSJeremy L Thompson 
365f80ce2aSJacob Faibussowitsch   PetscFunctionBeginUser;
375f80ce2aSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
389566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
399566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
407a17eebaSJed Brown   {
417a17eebaSJed Brown     IS              field_is;
427a17eebaSJed Brown     const PetscInt *fields;
437a17eebaSJed Brown     PetscInt        num_fields;
445f80ce2aSJacob Faibussowitsch 
459566063dSJacob Faibussowitsch     PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds));
46a2c9b50fSJeremy L Thompson     // Translate dm_field to ds_field
479566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(field_is, &fields));
489566063dSJacob Faibussowitsch     PetscCall(ISGetSize(field_is, &num_fields));
497a17eebaSJed Brown     for (PetscInt i=0; i<num_fields; i++) {
507a17eebaSJed Brown       if (dm_field == fields[i]) {
517a17eebaSJed Brown         ds_field = i;
52a2c9b50fSJeremy L Thompson         break;
53a2c9b50fSJeremy L Thompson       }
54a2c9b50fSJeremy L Thompson     }
559566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(field_is, &fields));
56a2c9b50fSJeremy L Thompson   }
5763a3b9bcSJacob Faibussowitsch   PetscCheck(ds_field != -1,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
58a2c9b50fSJeremy L Thompson 
59a2c9b50fSJeremy L Thompson   {
60a2c9b50fSJeremy L Thompson     PetscInt depth;
61a2c9b50fSJeremy L Thompson     DMLabel  depth_label;
62a2c9b50fSJeremy L Thompson     IS       depth_is;
635f80ce2aSJacob Faibussowitsch 
649566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepth(dm, &depth));
659566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
669566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is));
67a2c9b50fSJeremy L Thompson     if (domain_label) {
68a2c9b50fSJeremy L Thompson       IS domain_is;
695f80ce2aSJacob Faibussowitsch 
709566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is));
71a2c9b50fSJeremy L Thompson       if (domain_is) { // domainIS is non-empty
729566063dSJacob Faibussowitsch         PetscCall(ISIntersect(depth_is, domain_is, &iter_is));
739566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&domain_is));
74a2c9b50fSJeremy L Thompson       } else { // domainIS is NULL (empty)
75a2c9b50fSJeremy L Thompson         iter_is = NULL;
76a2c9b50fSJeremy L Thompson       }
779566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&depth_is));
78a2c9b50fSJeremy L Thompson     } else {
79a2c9b50fSJeremy L Thompson       iter_is = depth_is;
80a2c9b50fSJeremy L Thompson     }
81a2c9b50fSJeremy L Thompson     if (iter_is) {
829566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(iter_is, num_cells));
839566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(iter_is, &iter_indices));
84a2c9b50fSJeremy L Thompson     } else {
85a2c9b50fSJeremy L Thompson       *num_cells = 0;
86a2c9b50fSJeremy L Thompson       iter_indices = NULL;
87a2c9b50fSJeremy L Thompson     }
88a2c9b50fSJeremy L Thompson   }
89a2c9b50fSJeremy L Thompson 
90a2c9b50fSJeremy L Thompson   {
91a2c9b50fSJeremy L Thompson     PetscDualSpace dual_space;
92a2c9b50fSJeremy L Thompson     PetscInt       num_dual_basis_vectors;
935f80ce2aSJacob Faibussowitsch 
949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject*)&fe));
959566063dSJacob Faibussowitsch     PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
96*7c48043bSMatthew G. Knepley     PetscCheck(fe, PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " is invalid for DG coordinates", height);
979566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
989566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
999566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetNumComponents(dual_space, num_comp));
10063a3b9bcSJacob 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);
101a2c9b50fSJeremy L Thompson     *cell_size = num_dual_basis_vectors / *num_comp;
102a2c9b50fSJeremy L Thompson   }
103a2c9b50fSJeremy L Thompson   PetscInt restr_size = (*num_cells)*(*cell_size);
1049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(restr_size, &restr_indices));
105a2c9b50fSJeremy L Thompson   PetscInt cell_offset = 0;
106a2c9b50fSJeremy L Thompson 
107a2c9b50fSJeremy L Thompson   PetscInt P = (PetscInt) pow(*cell_size, 1.0 / (dim - height));
108a2c9b50fSJeremy L Thompson   for (PetscInt p = 0; p < *num_cells; p++) {
109a2c9b50fSJeremy L Thompson     PetscBool flip = PETSC_FALSE;
110a2c9b50fSJeremy L Thompson     PetscInt c = iter_indices[p];
111a2c9b50fSJeremy L Thompson     PetscInt num_indices, *indices;
112a2c9b50fSJeremy L Thompson     PetscInt field_offsets[17]; // max number of fields plus 1
1139566063dSJacob Faibussowitsch     PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
114a2c9b50fSJeremy L Thompson     if (height > 0) {
115a2c9b50fSJeremy L Thompson       PetscInt num_cells_support, num_faces, start = -1;
116a2c9b50fSJeremy L Thompson       const PetscInt *orients, *faces, *cells;
1179566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, c, &cells));
1189566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, c, &num_cells_support));
1195f80ce2aSJacob 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);
1209566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, cells[0], &faces));
1219566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, cells[0], &num_faces));
122a2c9b50fSJeremy L Thompson       for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
1235f80ce2aSJacob Faibussowitsch       PetscCheck(start >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c);
1249566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeOrientation(dm, cells[0], &orients));
125a2c9b50fSJeremy L Thompson       if (orients[start] < 0) flip = PETSC_TRUE;
126a2c9b50fSJeremy L Thompson     }
127a2c9b50fSJeremy L Thompson 
128a2c9b50fSJeremy L Thompson     for (PetscInt i = 0; i < *cell_size; i++) {
129a2c9b50fSJeremy L Thompson       PetscInt ii = i;
130a2c9b50fSJeremy L Thompson       if (flip) {
131a2c9b50fSJeremy L Thompson         if (*cell_size == P) ii = *cell_size - 1 - i;
132a2c9b50fSJeremy L Thompson         else if (*cell_size == P*P) {
133a2c9b50fSJeremy L Thompson           PetscInt row = i / P, col = i % P;
134a2c9b50fSJeremy L Thompson           ii = row + col * P;
13563a3b9bcSJacob 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);
136a2c9b50fSJeremy L Thompson       }
137a2c9b50fSJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
138a2c9b50fSJeremy L Thompson       PetscInt loc = indices[field_offsets[dm_field] + ii*(*num_comp)];
139a2c9b50fSJeremy L Thompson       restr_indices[cell_offset++] = loc >= 0 ? loc : -(loc + 1);
140a2c9b50fSJeremy L Thompson     }
1419566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL));
142a2c9b50fSJeremy L Thompson   }
14363a3b9bcSJacob 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);
1449566063dSJacob Faibussowitsch   if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices));
1459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iter_is));
146a2c9b50fSJeremy L Thompson 
147a2c9b50fSJeremy L Thompson   *offsets = restr_indices;
1489566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, l_size));
149a2c9b50fSJeremy L Thompson   PetscFunctionReturn(0);
150a2c9b50fSJeremy L Thompson }
151a2c9b50fSJeremy L Thompson 
152a2c9b50fSJeremy L Thompson #if defined(PETSC_HAVE_LIBCEED)
153f918ec44SMatthew G. Knepley #include <petscdmplexceed.h>
154f918ec44SMatthew G. Knepley 
155a2c9b50fSJeremy L Thompson /*@C
156a2c9b50fSJeremy L Thompson   DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector)
157a2c9b50fSJeremy L Thompson 
158a2c9b50fSJeremy L Thompson   Input Parameters:
159a2c9b50fSJeremy L Thompson   dm - The DMPlex object
160a2c9b50fSJeremy L Thompson   domain_label - label for DMPlex domain, or NULL for the whole domain
161a2c9b50fSJeremy L Thompson   label_value - Stratum value
162a2c9b50fSJeremy L Thompson   height - Height of target cells in DMPlex topology
163a2c9b50fSJeremy L Thompson   dm_field - Index of DMPlex field
164a2c9b50fSJeremy L Thompson 
165a2c9b50fSJeremy L Thompson   Output Parameters:
166a2c9b50fSJeremy L Thompson   ERestrict - libCEED restriction from local vector to to the cells
167a2c9b50fSJeremy L Thompson 
168a2c9b50fSJeremy L Thompson   Level: developer
169a2c9b50fSJeremy L Thompson @*/
170a2c9b50fSJeremy L Thompson PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict)
171f918ec44SMatthew G. Knepley {
172f918ec44SMatthew G. Knepley   PetscFunctionBeginUser;
173f918ec44SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
174f918ec44SMatthew G. Knepley   PetscValidPointer(ERestrict, 2);
175f918ec44SMatthew G. Knepley   if (!dm->ceedERestrict) {
176a2c9b50fSJeremy L Thompson     PetscInt     num_cells, cell_size, num_comp, lvec_size, *restr_indices;
177a2c9b50fSJeremy L Thompson     CeedElemRestriction elem_restr;
178f918ec44SMatthew G. Knepley     Ceed         ceed;
179f918ec44SMatthew G. Knepley 
1809566063dSJacob Faibussowitsch     PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices));
1819566063dSJacob Faibussowitsch     PetscCall(DMGetCeed(dm, &ceed));
1829566063dSJacob Faibussowitsch     PetscCallCEED(CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, restr_indices, &elem_restr));
1839566063dSJacob Faibussowitsch     PetscCall(PetscFree(restr_indices));
184a2c9b50fSJeremy L Thompson     dm->ceedERestrict = elem_restr;
185f918ec44SMatthew G. Knepley   }
186f918ec44SMatthew G. Knepley   *ERestrict = dm->ceedERestrict;
187f918ec44SMatthew G. Knepley   PetscFunctionReturn(0);
188f918ec44SMatthew G. Knepley }
189f918ec44SMatthew G. Knepley 
190f918ec44SMatthew G. Knepley #endif
191