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