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