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 PetscDS ds = NULL; 29a2c9b50fSJeremy L Thompson PetscFE fe; 30a2c9b50fSJeremy L Thompson PetscSection section; 317a17eebaSJed Brown PetscInt dim, ds_field = -1; 32a2c9b50fSJeremy L Thompson PetscInt *restr_indices; 33a2c9b50fSJeremy L Thompson const PetscInt *iter_indices; 34a2c9b50fSJeremy L Thompson IS iter_is; 35a2c9b50fSJeremy L Thompson 36*5f80ce2aSJacob Faibussowitsch PetscFunctionBeginUser; 37*5f80ce2aSJacob Faibussowitsch PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(dm, §ion)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 407a17eebaSJed Brown { 417a17eebaSJed Brown IS field_is; 427a17eebaSJed Brown const PetscInt *fields; 437a17eebaSJed Brown PetscInt num_fields; 44*5f80ce2aSJacob Faibussowitsch 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetRegionDS(dm, domain_label, &field_is, &ds)); 46a2c9b50fSJeremy L Thompson // Translate dm_field to ds_field 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(field_is, &fields)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetSize(field_is, &num_fields)); 497a17eebaSJed Brown for (PetscInt i=0; i<num_fields; i++) { 507a17eebaSJed Brown if (dm_field == fields[i]) { 517a17eebaSJed Brown ds_field = i; 52a2c9b50fSJeremy L Thompson break; 53a2c9b50fSJeremy L Thompson } 54a2c9b50fSJeremy L Thompson } 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(field_is, &fields)); 56a2c9b50fSJeremy L Thompson } 57*5f80ce2aSJacob Faibussowitsch PetscCheck(ds_field != -1,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Could not find dm_field %D in DS", dm_field); 58a2c9b50fSJeremy L Thompson 59a2c9b50fSJeremy L Thompson { 60a2c9b50fSJeremy L Thompson PetscInt depth; 61a2c9b50fSJeremy L Thompson DMLabel depth_label; 62a2c9b50fSJeremy L Thompson IS depth_is; 63*5f80ce2aSJacob Faibussowitsch 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepth(dm, &depth)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthLabel(dm, &depth_label)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 67a2c9b50fSJeremy L Thompson if (domain_label) { 68a2c9b50fSJeremy L Thompson IS domain_is; 69*5f80ce2aSJacob Faibussowitsch 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetStratumIS(domain_label, label_value, &domain_is)); 71a2c9b50fSJeremy L Thompson if (domain_is) { // domainIS is non-empty 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISIntersect(depth_is, domain_is, &iter_is)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&domain_is)); 74a2c9b50fSJeremy L Thompson } else { // domainIS is NULL (empty) 75a2c9b50fSJeremy L Thompson iter_is = NULL; 76a2c9b50fSJeremy L Thompson } 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&depth_is)); 78a2c9b50fSJeremy L Thompson } else { 79a2c9b50fSJeremy L Thompson iter_is = depth_is; 80a2c9b50fSJeremy L Thompson } 81a2c9b50fSJeremy L Thompson if (iter_is) { 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iter_is, num_cells)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iter_is, &iter_indices)); 84a2c9b50fSJeremy L Thompson } else { 85a2c9b50fSJeremy L Thompson *num_cells = 0; 86a2c9b50fSJeremy L Thompson iter_indices = NULL; 87a2c9b50fSJeremy L Thompson } 88a2c9b50fSJeremy L Thompson } 89a2c9b50fSJeremy L Thompson 90a2c9b50fSJeremy L Thompson { 91a2c9b50fSJeremy L Thompson PetscDualSpace dual_space; 92a2c9b50fSJeremy L Thompson PetscInt num_dual_basis_vectors; 93*5f80ce2aSJacob Faibussowitsch 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetDiscretization(ds, ds_field, (PetscObject*)&fe)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetHeightSubspace(fe, height, &fe)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetDualSpace(fe, &dual_space)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetNumComponents(dual_space, num_comp)); 99*5f80ce2aSJacob Faibussowitsch PetscCheck(num_dual_basis_vectors % *num_comp == 0,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); 100a2c9b50fSJeremy L Thompson *cell_size = num_dual_basis_vectors / *num_comp; 101a2c9b50fSJeremy L Thompson } 102a2c9b50fSJeremy L Thompson PetscInt restr_size = (*num_cells)*(*cell_size); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSupport(dm, c, &cells)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSupportSize(dm, c, &num_cells_support)); 118*5f80ce2aSJacob 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); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, cells[0], &faces)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeSize(dm, cells[0], &num_faces)); 121a2c9b50fSJeremy L Thompson for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;} 122*5f80ce2aSJacob Faibussowitsch PetscCheck(start >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %" PetscInt_FMT " in cone of its support", c); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeOrientation(dm, cells[0], &orients)); 124a2c9b50fSJeremy L Thompson if (orients[start] < 0) flip = PETSC_TRUE; 125a2c9b50fSJeremy L Thompson } 126a2c9b50fSJeremy L Thompson 127a2c9b50fSJeremy L Thompson for (PetscInt i = 0; i < *cell_size; i++) { 128a2c9b50fSJeremy L Thompson PetscInt ii = i; 129a2c9b50fSJeremy L Thompson if (flip) { 130a2c9b50fSJeremy L Thompson if (*cell_size == P) ii = *cell_size - 1 - i; 131a2c9b50fSJeremy L Thompson else if (*cell_size == P*P) { 132a2c9b50fSJeremy L Thompson PetscInt row = i / P, col = i % P; 133a2c9b50fSJeremy L Thompson ii = row + col * P; 13498921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for flipping point with cell size %D != P (%D) or P^2", *cell_size, P); 135a2c9b50fSJeremy L Thompson } 136a2c9b50fSJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 137a2c9b50fSJeremy L Thompson PetscInt loc = indices[field_offsets[dm_field] + ii*(*num_comp)]; 138a2c9b50fSJeremy L Thompson restr_indices[cell_offset++] = loc >= 0 ? loc : -(loc + 1); 139a2c9b50fSJeremy L Thompson } 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &num_indices, &indices, field_offsets, NULL)); 141a2c9b50fSJeremy L Thompson } 142*5f80ce2aSJacob Faibussowitsch PetscCheck(cell_offset == restr_size,PETSC_COMM_SELF, PETSC_ERR_SUP, "Shape mismatch, offsets array of shape (%D, %D) initialized for %D nodes", *num_cells, (*cell_size), cell_offset); 143*5f80ce2aSJacob Faibussowitsch if (iter_is) CHKERRQ(ISRestoreIndices(iter_is, &iter_indices)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iter_is)); 145a2c9b50fSJeremy L Thompson 146a2c9b50fSJeremy L Thompson *offsets = restr_indices; 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetStorageSize(section, l_size)); 148a2c9b50fSJeremy L Thompson PetscFunctionReturn(0); 149a2c9b50fSJeremy L Thompson } 150a2c9b50fSJeremy L Thompson 151a2c9b50fSJeremy L Thompson #if defined(PETSC_HAVE_LIBCEED) 152f918ec44SMatthew G. Knepley #include <petscdmplexceed.h> 153f918ec44SMatthew G. Knepley 154a2c9b50fSJeremy L Thompson /*@C 155a2c9b50fSJeremy L Thompson DMPlexGetCeedRestriction - Define the libCEED map from the local vector (Lvector) to the cells (Evector) 156a2c9b50fSJeremy L Thompson 157a2c9b50fSJeremy L Thompson Input Parameters: 158a2c9b50fSJeremy L Thompson dm - The DMPlex object 159a2c9b50fSJeremy L Thompson domain_label - label for DMPlex domain, or NULL for the whole domain 160a2c9b50fSJeremy L Thompson label_value - Stratum value 161a2c9b50fSJeremy L Thompson height - Height of target cells in DMPlex topology 162a2c9b50fSJeremy L Thompson dm_field - Index of DMPlex field 163a2c9b50fSJeremy L Thompson 164a2c9b50fSJeremy L Thompson Output Parameters: 165a2c9b50fSJeremy L Thompson ERestrict - libCEED restriction from local vector to to the cells 166a2c9b50fSJeremy L Thompson 167a2c9b50fSJeremy L Thompson Level: developer 168a2c9b50fSJeremy L Thompson @*/ 169a2c9b50fSJeremy L Thompson PetscErrorCode DMPlexGetCeedRestriction(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, CeedElemRestriction *ERestrict) 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 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_cells, &cell_size, &num_comp, &lvec_size, &restr_indices)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCeed(dm, &ceed)); 181*5f80ce2aSJacob Faibussowitsch CHKERRQ_CEED(CeedElemRestrictionCreate(ceed, num_cells, cell_size, num_comp, 1, lvec_size, CEED_MEM_HOST, CEED_COPY_VALUES, restr_indices, &elem_restr)); 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(restr_indices)); 183a2c9b50fSJeremy L Thompson dm->ceedERestrict = elem_restr; 184f918ec44SMatthew G. Knepley } 185f918ec44SMatthew G. Knepley *ERestrict = dm->ceedERestrict; 186f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 187f918ec44SMatthew G. Knepley } 188f918ec44SMatthew G. Knepley 189f918ec44SMatthew G. Knepley #endif 190