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: 7*a1cb98faSBarry Smith + dm - The `DMPLEX` object 8*a1cb98faSBarry Smith . domain_label - label for `DMPLEX` domain, or NULL for whole domain 9*a1cb98faSBarry Smith . label_value - Stratum value 10*a1cb98faSBarry Smith . height - Height of target cells in `DMPLEX` topology 11*a1cb98faSBarry Smith - dm_field - Index of `DMPLEX` field 12a2c9b50fSJeremy L Thompson 13a2c9b50fSJeremy L Thompson Output Parameters: 14*a1cb98faSBarry Smith + num_cells - Number of local cells 15*a1cb98faSBarry Smith . cell_size - Size of each cell, given by cell_size * num_comp = num_dof 16*a1cb98faSBarry Smith . num_comp - Number of components per dof 17*a1cb98faSBarry Smith . l_size - Size of local vector 18*a1cb98faSBarry Smith - offsets - Allocated offsets array for cells 19a2c9b50fSJeremy L Thompson 20a2c9b50fSJeremy L Thompson Level: developer 21a2c9b50fSJeremy L Thompson 22*a1cb98faSBarry Smith Notes: 23*a1cb98faSBarry Smith 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]. 24*a1cb98faSBarry Smith 25*a1cb98faSBarry Smith Caller is responsible for freeing the offsets array using `PetscFree()`. 26*a1cb98faSBarry Smith 27*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetClosureIndices()`, `DMPlexSetClosurePermutationTensor()`, `DMPlexGetCeedRestriction()` 28a2c9b50fSJeremy L Thompson @*/ 29d71ae5a4SJacob Faibussowitsch 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) 30d71ae5a4SJacob Faibussowitsch { 31a2c9b50fSJeremy L Thompson PetscDS ds = NULL; 32a2c9b50fSJeremy L Thompson PetscFE fe; 33a2c9b50fSJeremy L Thompson PetscSection section; 347a17eebaSJed Brown PetscInt dim, ds_field = -1; 35a2c9b50fSJeremy L Thompson PetscInt *restr_indices; 36a2c9b50fSJeremy L Thompson const PetscInt *iter_indices; 37a2c9b50fSJeremy L Thompson IS iter_is; 38a2c9b50fSJeremy L Thompson 395f80ce2aSJacob Faibussowitsch PetscFunctionBeginUser; 405f80ce2aSJacob Faibussowitsch PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 419566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 429566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 437a17eebaSJed Brown { 447a17eebaSJed Brown IS field_is; 457a17eebaSJed Brown const PetscInt *fields; 467a17eebaSJed Brown PetscInt num_fields; 475f80ce2aSJacob Faibussowitsch 489566063dSJacob Faibussowitsch PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds)); 49a2c9b50fSJeremy L Thompson // Translate dm_field to ds_field 509566063dSJacob Faibussowitsch PetscCall(ISGetIndices(field_is, &fields)); 519566063dSJacob Faibussowitsch PetscCall(ISGetSize(field_is, &num_fields)); 527a17eebaSJed Brown for (PetscInt i = 0; i < num_fields; i++) { 537a17eebaSJed Brown if (dm_field == fields[i]) { 547a17eebaSJed Brown ds_field = i; 55a2c9b50fSJeremy L Thompson break; 56a2c9b50fSJeremy L Thompson } 57a2c9b50fSJeremy L Thompson } 589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(field_is, &fields)); 59a2c9b50fSJeremy L Thompson } 6063a3b9bcSJacob Faibussowitsch PetscCheck(ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 61a2c9b50fSJeremy L Thompson 62a2c9b50fSJeremy L Thompson { 63a2c9b50fSJeremy L Thompson PetscInt depth; 64a2c9b50fSJeremy L Thompson DMLabel depth_label; 65a2c9b50fSJeremy L Thompson IS depth_is; 665f80ce2aSJacob Faibussowitsch 679566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 689566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 699566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 70a2c9b50fSJeremy L Thompson if (domain_label) { 71a2c9b50fSJeremy L Thompson IS domain_is; 725f80ce2aSJacob Faibussowitsch 739566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is)); 74a2c9b50fSJeremy L Thompson if (domain_is) { // domainIS is non-empty 759566063dSJacob Faibussowitsch PetscCall(ISIntersect(depth_is, domain_is, &iter_is)); 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&domain_is)); 77a2c9b50fSJeremy L Thompson } else { // domainIS is NULL (empty) 78a2c9b50fSJeremy L Thompson iter_is = NULL; 79a2c9b50fSJeremy L Thompson } 809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&depth_is)); 81a2c9b50fSJeremy L Thompson } else { 82a2c9b50fSJeremy L Thompson iter_is = depth_is; 83a2c9b50fSJeremy L Thompson } 84a2c9b50fSJeremy L Thompson if (iter_is) { 859566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iter_is, num_cells)); 869566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iter_is, &iter_indices)); 87a2c9b50fSJeremy L Thompson } else { 88a2c9b50fSJeremy L Thompson *num_cells = 0; 89a2c9b50fSJeremy L Thompson iter_indices = NULL; 90a2c9b50fSJeremy L Thompson } 91a2c9b50fSJeremy L Thompson } 92a2c9b50fSJeremy L Thompson 93a2c9b50fSJeremy L Thompson { 94a2c9b50fSJeremy L Thompson PetscDualSpace dual_space; 95a2c9b50fSJeremy L Thompson PetscInt num_dual_basis_vectors; 965f80ce2aSJacob Faibussowitsch 979566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 989566063dSJacob Faibussowitsch PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 997c48043bSMatthew G. Knepley PetscCheck(fe, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " is invalid for DG coordinates", height); 1009566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 1019566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 1029566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(dual_space, num_comp)); 10363a3b9bcSJacob 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); 104a2c9b50fSJeremy L Thompson *cell_size = num_dual_basis_vectors / *num_comp; 105a2c9b50fSJeremy L Thompson } 106a2c9b50fSJeremy L Thompson PetscInt restr_size = (*num_cells) * (*cell_size); 1079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(restr_size, &restr_indices)); 108a2c9b50fSJeremy L Thompson PetscInt cell_offset = 0; 109a2c9b50fSJeremy L Thompson 110a2c9b50fSJeremy L Thompson PetscInt P = (PetscInt)pow(*cell_size, 1.0 / (dim - height)); 111a2c9b50fSJeremy L Thompson for (PetscInt p = 0; p < *num_cells; p++) { 112a2c9b50fSJeremy L Thompson PetscBool flip = PETSC_FALSE; 113a2c9b50fSJeremy L Thompson PetscInt c = iter_indices[p]; 114a2c9b50fSJeremy L Thompson PetscInt num_indices, *indices; 115a2c9b50fSJeremy L Thompson PetscInt field_offsets[17]; // max number of fields plus 1 1169566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL)); 117a2c9b50fSJeremy L Thompson if (height > 0) { 118a2c9b50fSJeremy L Thompson PetscInt num_cells_support, num_faces, start = -1; 119a2c9b50fSJeremy L Thompson const PetscInt *orients, *faces, *cells; 1209566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, c, &cells)); 1219566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, c, &num_cells_support)); 1225f80ce2aSJacob 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); 1239566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cells[0], &faces)); 1249566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &num_faces)); 1259371c9d4SSatish Balay for (PetscInt i = 0; i < num_faces; i++) { 1269371c9d4SSatish Balay if (faces[i] == c) start = i; 1279371c9d4SSatish Balay } 1285f80ce2aSJacob Faibussowitsch PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c); 1299566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, cells[0], &orients)); 130a2c9b50fSJeremy L Thompson if (orients[start] < 0) flip = PETSC_TRUE; 131a2c9b50fSJeremy L Thompson } 132a2c9b50fSJeremy L Thompson 133a2c9b50fSJeremy L Thompson for (PetscInt i = 0; i < *cell_size; i++) { 134a2c9b50fSJeremy L Thompson PetscInt ii = i; 135a2c9b50fSJeremy L Thompson if (flip) { 136a2c9b50fSJeremy L Thompson if (*cell_size == P) ii = *cell_size - 1 - i; 137a2c9b50fSJeremy L Thompson else if (*cell_size == P * P) { 138a2c9b50fSJeremy L Thompson PetscInt row = i / P, col = i % P; 139a2c9b50fSJeremy L Thompson ii = row + col * P; 14063a3b9bcSJacob 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); 141a2c9b50fSJeremy L Thompson } 142a2c9b50fSJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 143a2c9b50fSJeremy L Thompson PetscInt loc = indices[field_offsets[dm_field] + ii * (*num_comp)]; 144a2c9b50fSJeremy L Thompson restr_indices[cell_offset++] = loc >= 0 ? loc : -(loc + 1); 145a2c9b50fSJeremy L Thompson } 1469566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL)); 147a2c9b50fSJeremy L Thompson } 14863a3b9bcSJacob 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); 1499566063dSJacob Faibussowitsch if (iter_is) PetscCall(ISRestoreIndices(iter_is, &iter_indices)); 1509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iter_is)); 151a2c9b50fSJeremy L Thompson 152a2c9b50fSJeremy L Thompson *offsets = restr_indices; 1539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, l_size)); 154a2c9b50fSJeremy L Thompson PetscFunctionReturn(0); 155a2c9b50fSJeremy L Thompson } 156a2c9b50fSJeremy L Thompson 157a2c9b50fSJeremy L Thompson #if defined(PETSC_HAVE_LIBCEED) 158f918ec44SMatthew G. Knepley #include <petscdmplexceed.h> 159f918ec44SMatthew G. Knepley 160a2c9b50fSJeremy L Thompson /*@C 161a2c9b50fSJeremy L Thompson DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector) 162a2c9b50fSJeremy L Thompson 163a2c9b50fSJeremy L Thompson Input Parameters: 164*a1cb98faSBarry Smith + dm - The `DMPLEX` object 165*a1cb98faSBarry Smith . domain_label - label for `DMPLEX` domain, or NULL for the whole domain 166*a1cb98faSBarry Smith . label_value - Stratum value 167*a1cb98faSBarry Smith . height - Height of target cells in `DMPLEX` topology 168*a1cb98faSBarry Smith - dm_field - Index of `DMPLEX` field 169a2c9b50fSJeremy L Thompson 170*a1cb98faSBarry Smith Output Parameter: 171*a1cb98faSBarry Smith . ERestrict - libCEED restriction from local vector to to the cells 172a2c9b50fSJeremy L Thompson 173a2c9b50fSJeremy L Thompson Level: developer 174*a1cb98faSBarry Smith 175*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexGetLocalOffsets()` 176a2c9b50fSJeremy L Thompson @*/ 177d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict) 178d71ae5a4SJacob Faibussowitsch { 179f918ec44SMatthew G. Knepley PetscFunctionBeginUser; 180f918ec44SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 181f918ec44SMatthew G. Knepley PetscValidPointer(ERestrict, 2); 182f918ec44SMatthew G. Knepley if (!dm->ceedERestrict) { 183a2c9b50fSJeremy L Thompson PetscInt num_cells, cell_size, num_comp, lvec_size, *restr_indices; 184a2c9b50fSJeremy L Thompson CeedElemRestriction elem_restr; 185f918ec44SMatthew G. Knepley Ceed ceed; 186f918ec44SMatthew G. Knepley 1879566063dSJacob Faibussowitsch PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetCeed(dm, &ceed)); 1899566063dSJacob Faibussowitsch PetscCallCEED(CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, restr_indices, &elem_restr)); 1909566063dSJacob Faibussowitsch PetscCall(PetscFree(restr_indices)); 191a2c9b50fSJeremy L Thompson dm->ceedERestrict = elem_restr; 192f918ec44SMatthew G. Knepley } 193f918ec44SMatthew G. Knepley *ERestrict = dm->ceedERestrict; 194f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 195f918ec44SMatthew G. Knepley } 196f918ec44SMatthew G. Knepley 197f918ec44SMatthew G. Knepley #endif 198