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