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 @*/ 26*9371c9d4SSatish Balay 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 PetscDS ds = NULL; 28a2c9b50fSJeremy L Thompson PetscFE fe; 29a2c9b50fSJeremy L Thompson PetscSection section; 307a17eebaSJed Brown PetscInt dim, ds_field = -1; 31a2c9b50fSJeremy L Thompson PetscInt *restr_indices; 32a2c9b50fSJeremy L Thompson const PetscInt *iter_indices; 33a2c9b50fSJeremy L Thompson IS iter_is; 34a2c9b50fSJeremy L Thompson 355f80ce2aSJacob Faibussowitsch PetscFunctionBeginUser; 365f80ce2aSJacob Faibussowitsch PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 379566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 389566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 397a17eebaSJed Brown { 407a17eebaSJed Brown IS field_is; 417a17eebaSJed Brown const PetscInt *fields; 427a17eebaSJed Brown PetscInt num_fields; 435f80ce2aSJacob Faibussowitsch 449566063dSJacob Faibussowitsch PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds)); 45a2c9b50fSJeremy L Thompson // Translate dm_field to ds_field 469566063dSJacob Faibussowitsch PetscCall(ISGetIndices(field_is, &fields)); 479566063dSJacob Faibussowitsch PetscCall(ISGetSize(field_is, &num_fields)); 487a17eebaSJed Brown for (PetscInt i = 0; i < num_fields; i++) { 497a17eebaSJed Brown if (dm_field == fields[i]) { 507a17eebaSJed Brown ds_field = i; 51a2c9b50fSJeremy L Thompson break; 52a2c9b50fSJeremy L Thompson } 53a2c9b50fSJeremy L Thompson } 549566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(field_is, &fields)); 55a2c9b50fSJeremy L Thompson } 5663a3b9bcSJacob Faibussowitsch PetscCheck(ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 57a2c9b50fSJeremy L Thompson 58a2c9b50fSJeremy L Thompson { 59a2c9b50fSJeremy L Thompson PetscInt depth; 60a2c9b50fSJeremy L Thompson DMLabel depth_label; 61a2c9b50fSJeremy L Thompson IS depth_is; 625f80ce2aSJacob Faibussowitsch 639566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 649566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 659566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 66a2c9b50fSJeremy L Thompson if (domain_label) { 67a2c9b50fSJeremy L Thompson IS domain_is; 685f80ce2aSJacob Faibussowitsch 699566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is)); 70a2c9b50fSJeremy L Thompson if (domain_is) { // domainIS is non-empty 719566063dSJacob Faibussowitsch PetscCall(ISIntersect(depth_is, domain_is, &iter_is)); 729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&domain_is)); 73a2c9b50fSJeremy L Thompson } else { // domainIS is NULL (empty) 74a2c9b50fSJeremy L Thompson iter_is = NULL; 75a2c9b50fSJeremy L Thompson } 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&depth_is)); 77a2c9b50fSJeremy L Thompson } else { 78a2c9b50fSJeremy L Thompson iter_is = depth_is; 79a2c9b50fSJeremy L Thompson } 80a2c9b50fSJeremy L Thompson if (iter_is) { 819566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iter_is, num_cells)); 829566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iter_is, &iter_indices)); 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; 925f80ce2aSJacob Faibussowitsch 939566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 949566063dSJacob Faibussowitsch PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 957c48043bSMatthew G. Knepley PetscCheck(fe, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " is invalid for DG coordinates", height); 969566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 979566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 989566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(dual_space, num_comp)); 9963a3b9bcSJacob 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); 100a2c9b50fSJeremy L Thompson *cell_size = num_dual_basis_vectors / *num_comp; 101a2c9b50fSJeremy L Thompson } 102a2c9b50fSJeremy L Thompson PetscInt restr_size = (*num_cells) * (*cell_size); 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(restr_size, &restr_indices)); 104a2c9b50fSJeremy L Thompson PetscInt cell_offset = 0; 105a2c9b50fSJeremy L Thompson 106a2c9b50fSJeremy L Thompson PetscInt P = (PetscInt)pow(*cell_size, 1.0 / (dim - height)); 107a2c9b50fSJeremy L Thompson for (PetscInt p = 0; p < *num_cells; p++) { 108a2c9b50fSJeremy L Thompson PetscBool flip = PETSC_FALSE; 109a2c9b50fSJeremy L Thompson PetscInt c = iter_indices[p]; 110a2c9b50fSJeremy L Thompson PetscInt num_indices, *indices; 111a2c9b50fSJeremy L Thompson PetscInt field_offsets[17]; // max number of fields plus 1 1129566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL)); 113a2c9b50fSJeremy L Thompson if (height > 0) { 114a2c9b50fSJeremy L Thompson PetscInt num_cells_support, num_faces, start = -1; 115a2c9b50fSJeremy L Thompson const PetscInt *orients, *faces, *cells; 1169566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, c, &cells)); 1179566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, c, &num_cells_support)); 1185f80ce2aSJacob 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); 1199566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cells[0], &faces)); 1209566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &num_faces)); 121*9371c9d4SSatish Balay for (PetscInt i = 0; i < num_faces; i++) { 122*9371c9d4SSatish Balay if (faces[i] == c) start = i; 123*9371c9d4SSatish Balay } 1245f80ce2aSJacob Faibussowitsch PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c); 1259566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, cells[0], &orients)); 126a2c9b50fSJeremy L Thompson if (orients[start] < 0) flip = PETSC_TRUE; 127a2c9b50fSJeremy L Thompson } 128a2c9b50fSJeremy L Thompson 129a2c9b50fSJeremy L Thompson for (PetscInt i = 0; i < *cell_size; i++) { 130a2c9b50fSJeremy L Thompson PetscInt ii = i; 131a2c9b50fSJeremy L Thompson if (flip) { 132a2c9b50fSJeremy L Thompson if (*cell_size == P) ii = *cell_size - 1 - i; 133a2c9b50fSJeremy L Thompson else if (*cell_size == P * P) { 134a2c9b50fSJeremy L Thompson PetscInt row = i / P, col = i % P; 135a2c9b50fSJeremy L Thompson ii = row + col * P; 13663a3b9bcSJacob 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); 137a2c9b50fSJeremy L Thompson } 138a2c9b50fSJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 139a2c9b50fSJeremy L Thompson PetscInt loc = indices[field_offsets[dm_field] + ii * (*num_comp)]; 140a2c9b50fSJeremy L Thompson restr_indices[cell_offset++] = loc >= 0 ? loc : -(loc + 1); 141a2c9b50fSJeremy L Thompson } 1429566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL)); 143a2c9b50fSJeremy L Thompson } 14463a3b9bcSJacob 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); 1459566063dSJacob Faibussowitsch if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices)); 1469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iter_is)); 147a2c9b50fSJeremy L Thompson 148a2c9b50fSJeremy L Thompson *offsets = restr_indices; 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, l_size)); 150a2c9b50fSJeremy L Thompson PetscFunctionReturn(0); 151a2c9b50fSJeremy L Thompson } 152a2c9b50fSJeremy L Thompson 153a2c9b50fSJeremy L Thompson #if defined(PETSC_HAVE_LIBCEED) 154f918ec44SMatthew G. Knepley #include <petscdmplexceed.h> 155f918ec44SMatthew G. Knepley 156a2c9b50fSJeremy L Thompson /*@C 157a2c9b50fSJeremy L Thompson DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector) 158a2c9b50fSJeremy L Thompson 159a2c9b50fSJeremy L Thompson Input Parameters: 160a2c9b50fSJeremy L Thompson dm - The DMPlex object 161a2c9b50fSJeremy L Thompson domain_label - label for DMPlex domain, or NULL for the whole domain 162a2c9b50fSJeremy L Thompson label_value - Stratum value 163a2c9b50fSJeremy L Thompson height - Height of target cells in DMPlex topology 164a2c9b50fSJeremy L Thompson dm_field - Index of DMPlex field 165a2c9b50fSJeremy L Thompson 166a2c9b50fSJeremy L Thompson Output Parameters: 167a2c9b50fSJeremy L Thompson ERestrict - libCEED restriction from local vector to to the cells 168a2c9b50fSJeremy L Thompson 169a2c9b50fSJeremy L Thompson Level: developer 170a2c9b50fSJeremy L Thompson @*/ 171*9371c9d4SSatish Balay PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict) { 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