xref: /petsc/src/dm/impls/plex/plexceed.c (revision db2b9530cea7203cdddcdfcaf19eada576ef3459)
1 #include <petsc/private/dmpleximpl.h>           /*I      "petscdmplex.h"          I*/
2 
3 /*@C
4   DMPlexGetLocalOffsets - Allocate and populate array of local offsets.
5 
6   Input Parameters:
7   dm - The DMPlex object
8   domain_label - label for DMPlex domain, or NULL for whole domain
9   label_value - Stratum value
10   height - Height of target cells in DMPlex topology
11   dm_field - Index of DMPlex field
12 
13   Output Parameters:
14   num_cells - Number of local cells
15   cell_size - Size of each cell, given by cell_size * num_comp = num_dof
16   num_comp - Number of components per dof
17   l_size - Size of local vector
18   offsets - Allocated offsets array for cells
19 
20   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().
21 
22   Level: developer
23 
24 .seealso: DMPlexGetClosureIndices(), DMPlexSetClosurePermutationTensor(), DMPlexGetCeedRestriction()
25 @*/
26 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)
27 {
28   PetscErrorCode ierr;
29 
30   PetscFunctionBeginUser;
31   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
32   PetscDS      ds = NULL;
33   PetscFE      fe;
34   PetscSection section;
35   PetscInt     dim, ds_field = -1;
36   PetscInt    *restr_indices;
37   const PetscInt *iter_indices;
38   IS           iter_is;
39 
40   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
41   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
42   {
43     IS field_is;
44     const PetscInt *fields;
45     PetscInt num_fields;
46     ierr = DMGetRegionDS(dm, domain_label, &field_is, &ds);CHKERRQ(ierr);
47     // Translate dm_field to ds_field
48     ierr = ISGetIndices(field_is, &fields);CHKERRQ(ierr);
49     ierr = ISGetSize(field_is, &num_fields);CHKERRQ(ierr);
50     for (PetscInt i=0; i<num_fields; i++) {
51       if (dm_field == fields[i]) {
52         ds_field = i;
53         break;
54       }
55     }
56     ierr = ISRestoreIndices(field_is, &fields);CHKERRQ(ierr);
57   }
58   PetscCheckFalse(ds_field == -1,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Could not find dm_field %D in DS", dm_field);
59 
60   {
61     PetscInt depth;
62     DMLabel depth_label;
63     IS depth_is;
64     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
65     ierr = DMPlexGetDepthLabel(dm, &depth_label);CHKERRQ(ierr);
66     ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is);CHKERRQ(ierr);
67     if (domain_label) {
68       IS domain_is;
69       ierr = DMLabelGetStratumIS(domain_label, label_value, &domain_is);CHKERRQ(ierr);
70       if (domain_is) { // domainIS is non-empty
71         ierr = ISIntersect(depth_is, domain_is, &iter_is);CHKERRQ(ierr);
72         ierr = ISDestroy(&domain_is);CHKERRQ(ierr);
73       } else { // domainIS is NULL (empty)
74         iter_is = NULL;
75       }
76       ierr = ISDestroy(&depth_is);CHKERRQ(ierr);
77     } else {
78       iter_is = depth_is;
79     }
80     if (iter_is) {
81       ierr = ISGetLocalSize(iter_is, num_cells);CHKERRQ(ierr);
82       ierr = ISGetIndices(iter_is, &iter_indices);CHKERRQ(ierr);
83     } else {
84       *num_cells = 0;
85       iter_indices = NULL;
86     }
87   }
88 
89   {
90     PetscDualSpace dual_space;
91     PetscInt num_dual_basis_vectors;
92     ierr = PetscDSGetDiscretization(ds, ds_field, (PetscObject*)&fe);CHKERRQ(ierr);
93     ierr = PetscFEGetHeightSubspace(fe, height, &fe);CHKERRQ(ierr);
94     ierr = PetscFEGetDualSpace(fe, &dual_space);CHKERRQ(ierr);
95     ierr = PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors);CHKERRQ(ierr);
96     ierr = PetscDualSpaceGetNumComponents(dual_space, num_comp);CHKERRQ(ierr);
97     PetscCheckFalse(num_dual_basis_vectors % *num_comp != 0,PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for number of dual basis vectors %D not divisible by %D components", num_dual_basis_vectors, *num_comp);
98     *cell_size = num_dual_basis_vectors / *num_comp;
99   }
100   PetscInt restr_size = (*num_cells)*(*cell_size);
101   ierr = PetscMalloc1(restr_size, &restr_indices);CHKERRQ(ierr);
102   PetscInt cell_offset = 0;
103 
104   PetscInt P = (PetscInt) pow(*cell_size, 1.0 / (dim - height));
105   for (PetscInt p = 0; p < *num_cells; p++) {
106     PetscBool flip = PETSC_FALSE;
107     PetscInt c = iter_indices[p];
108     PetscInt num_indices, *indices;
109     PetscInt field_offsets[17]; // max number of fields plus 1
110     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL);CHKERRQ(ierr);
111     if (height > 0) {
112       PetscInt num_cells_support, num_faces, start = -1;
113       const PetscInt *orients, *faces, *cells;
114       ierr = DMPlexGetSupport(dm, c, &cells);CHKERRQ(ierr);
115       ierr = DMPlexGetSupportSize(dm, c, &num_cells_support);CHKERRQ(ierr);
116       PetscCheckFalse(num_cells_support != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Expected one cell in support of exterior face, but got %D cells", num_cells_support);
117       ierr = DMPlexGetCone(dm, cells[0], &faces);CHKERRQ(ierr);
118       ierr = DMPlexGetConeSize(dm, cells[0], &num_faces);CHKERRQ(ierr);
119       for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
120       PetscCheckFalse(start < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %D in cone of its support", c);
121       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients);CHKERRQ(ierr);
122       if (orients[start] < 0) flip = PETSC_TRUE;
123     }
124 
125     for (PetscInt i = 0; i < *cell_size; i++) {
126       PetscInt ii = i;
127       if (flip) {
128         if (*cell_size == P) ii = *cell_size - 1 - i;
129         else if (*cell_size == P*P) {
130           PetscInt row = i / P, col = i % P;
131           ii = row + col * P;
132         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for flipping point with cell size %D != P (%D) or P^2", *cell_size, P);
133       }
134       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
135       PetscInt loc = indices[field_offsets[dm_field] + ii*(*num_comp)];
136       restr_indices[cell_offset++] = loc >= 0 ? loc : -(loc + 1);
137     }
138     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL);CHKERRQ(ierr);
139   }
140   PetscCheckFalse(cell_offset != restr_size,PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, offsets array of shape (%D, %D) initialized for %D nodes", *num_cells, (*cell_size), cell_offset);
141   if (iter_is) { ierr = ISRestoreIndices(iter_is, &iter_indices);CHKERRQ(ierr); }
142   ierr = ISDestroy(&iter_is);CHKERRQ(ierr);
143 
144   *offsets = restr_indices;
145   ierr = PetscSectionGetStorageSize(section, l_size);CHKERRQ(ierr);
146   PetscFunctionReturn(0);
147 }
148 
149 #if defined(PETSC_HAVE_LIBCEED)
150 #include <petscdmplexceed.h>
151 
152 /*@C
153   DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector)
154 
155   Input Parameters:
156   dm - The DMPlex object
157   domain_label - label for DMPlex domain, or NULL for the whole domain
158   label_value - Stratum value
159   height - Height of target cells in DMPlex topology
160   dm_field - Index of DMPlex field
161 
162   Output Parameters:
163   ERestrict - libCEED restriction from local vector to to the cells
164 
165   Level: developer
166 @*/
167 PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict)
168 {
169   PetscErrorCode ierr;
170 
171   PetscFunctionBeginUser;
172   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
173   PetscValidPointer(ERestrict, 2);
174   if (!dm->ceedERestrict) {
175     PetscInt     num_cells, cell_size, num_comp, lvec_size, *restr_indices;
176     CeedElemRestriction elem_restr;
177     Ceed         ceed;
178 
179     ierr = DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices);CHKERRQ(ierr);
180     ierr = DMGetCeed(dm, &ceed);CHKERRQ(ierr);
181     ierr = CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, restr_indices, &elem_restr);CHKERRQ_CEED(ierr);
182     ierr = PetscFree(restr_indices);CHKERRQ(ierr);
183     dm->ceedERestrict = elem_restr;
184   }
185   *ERestrict = dm->ceedERestrict;
186   PetscFunctionReturn(0);
187 }
188 
189 #endif
190