1 #include "../include/petscutils.h" 2 #include "../include/petscmacros.h" 3 4 // ----------------------------------------------------------------------------- 5 // Convert PETSc MemType to libCEED MemType 6 // ----------------------------------------------------------------------------- 7 CeedMemType MemTypeP2C(PetscMemType mem_type) { 8 return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; 9 } 10 11 // ----------------------------------------------------------------------------- 12 // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c 13 // ----------------------------------------------------------------------------- 14 PetscErrorCode ProjectToUnitSphere(DM dm) { 15 Vec coordinates; 16 PetscScalar *coords; 17 PetscInt Nv, v, dim, d; 18 PetscErrorCode ierr; 19 20 PetscFunctionBeginUser; 21 ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr); 22 ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr); 23 ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr); 24 Nv /= dim; 25 ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr); 26 for (v = 0; v < Nv; ++v) { 27 PetscReal r = 0.0; 28 29 for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d])); 30 r = PetscSqrtReal(r); 31 for (d = 0; d < dim; ++d) coords[v*dim+d] /= r; 32 } 33 ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr); 34 PetscFunctionReturn(0); 35 }; 36 37 // ----------------------------------------------------------------------------- 38 // Apply 3D Kershaw mesh transformation 39 // ----------------------------------------------------------------------------- 40 // Transition from a value of "a" for x=0, to a value of "b" for x=1. Optionally 41 // smooth -- see the commented versions at the end. 42 static double step(const double a, const double b, double x) { 43 if (x <= 0) return a; 44 if (x >= 1) return b; 45 return a + (b-a) * (x); 46 } 47 48 // 1D transformation at the right boundary 49 static double right(const double eps, const double x) { 50 return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1); 51 } 52 53 // 1D transformation at the left boundary 54 static double left(const double eps, const double x) { 55 return 1-right(eps,1-x); 56 } 57 58 // Apply 3D Kershaw mesh transformation 59 // The eps parameters are in (0, 1] 60 // Uniform mesh is recovered for eps=1 61 PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) { 62 PetscErrorCode ierr; 63 Vec coord; 64 PetscInt ncoord; 65 PetscScalar *c; 66 67 PetscFunctionBeginUser; 68 ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr); 69 ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr); 70 ierr = VecGetArray(coord, &c); CHKERRQ(ierr); 71 72 for (PetscInt i = 0; i < ncoord; i += 3) { 73 PetscScalar x = c[i], y = c[i+1], z = c[i+2]; 74 PetscInt layer = x*6; 75 PetscScalar lambda = (x-layer/6.0)*6; 76 c[i] = x; 77 78 switch (layer) { 79 case 0: 80 c[i+1] = left(eps, y); 81 c[i+2] = left(eps, z); 82 break; 83 case 1: 84 case 4: 85 c[i+1] = step(left(eps, y), right(eps, y), lambda); 86 c[i+2] = step(left(eps, z), right(eps, z), lambda); 87 break; 88 case 2: 89 c[i+1] = step(right(eps, y), left(eps, y), lambda/2); 90 c[i+2] = step(right(eps, z), left(eps, z), lambda/2); 91 break; 92 case 3: 93 c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2); 94 c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2); 95 break; 96 default: 97 c[i+1] = right(eps, y); 98 c[i+2] = right(eps, z); 99 } 100 } 101 ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr); 102 PetscFunctionReturn(0); 103 } 104 105 // ----------------------------------------------------------------------------- 106 // Create BC label 107 // ----------------------------------------------------------------------------- 108 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 109 int ierr; 110 DMLabel label; 111 112 PetscFunctionBeginUser; 113 114 ierr = DMCreateLabel(dm, name); CHKERRQ(ierr); 115 ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr); 116 ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr); 117 ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr); 118 119 PetscFunctionReturn(0); 120 }; 121 122 // ----------------------------------------------------------------------------- 123 // This function sets up a DM for a given degree 124 // ----------------------------------------------------------------------------- 125 PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u, 126 PetscInt dim, bool enforce_bc, BCFunction bc_func) { 127 PetscInt ierr, marker_ids[1] = {1}; 128 PetscFE fe; 129 MPI_Comm comm; 130 131 PetscFunctionBeginUser; 132 133 // Setup FE 134 ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr); 135 ierr = PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, degree, degree, 136 &fe); CHKERRQ(ierr); 137 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 138 ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 139 140 // Setup DM 141 ierr = DMCreateDS(dm); CHKERRQ(ierr); 142 if (enforce_bc) { 143 PetscBool has_label; 144 DMHasLabel(dm, "marker", &has_label); 145 if (!has_label) {CreateBCLabel(dm, "marker");} 146 DMLabel label; 147 ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr); 148 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, "marker", 1, 149 marker_ids, 0, 0, NULL, 150 (void(*)(void))bc_func, NULL, NULL, NULL); 151 CHKERRQ(ierr); 152 } 153 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 154 CHKERRQ(ierr); 155 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 156 157 PetscFunctionReturn(0); 158 }; 159 160 // ----------------------------------------------------------------------------- 161 // Utility function - essential BC dofs are encoded in closure indices as -(i+1) 162 // ----------------------------------------------------------------------------- 163 PetscInt Involute(PetscInt i) { 164 return i >= 0 ? i : -(i + 1); 165 }; 166 167 // ----------------------------------------------------------------------------- 168 // Get CEED restriction data from DMPlex 169 // ----------------------------------------------------------------------------- 170 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 171 CeedInt topo_dim, CeedInt height, DMLabel domain_label, CeedInt value, 172 CeedElemRestriction *elem_restr) { 173 PetscSection section; 174 PetscInt p, num_elem, num_dof, *elem_restr_offsets, e_offset, num_fields, dim, 175 depth; 176 DMLabel depth_label; 177 IS depth_is, iter_is; 178 Vec U_loc; 179 const PetscInt *iter_indices; 180 PetscErrorCode ierr; 181 182 PetscFunctionBeginUser; 183 184 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 185 dim -= height; 186 ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 187 ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr); 188 PetscInt num_comp[num_fields], field_off[num_fields+1]; 189 field_off[0] = 0; 190 for (PetscInt f = 0; f < num_fields; f++) { 191 ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr); 192 field_off[f+1] = field_off[f] + num_comp[f]; 193 } 194 195 ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 196 ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr); 197 ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is); 198 CHKERRQ(ierr); 199 if (domain_label) { 200 IS domain_is; 201 ierr = DMLabelGetStratumIS(domain_label, value, &domain_is); CHKERRQ(ierr); 202 if (domain_is) { // domain_is is non-empty 203 ierr = ISIntersect(depth_is, domain_is, &iter_is); CHKERRQ(ierr); 204 ierr = ISDestroy(&domain_is); CHKERRQ(ierr); 205 } else { // domain_is is NULL (empty) 206 iter_is = NULL; 207 } 208 ierr = ISDestroy(&depth_is); CHKERRQ(ierr); 209 } else { 210 iter_is = depth_is; 211 } 212 if (iter_is) { 213 ierr = ISGetLocalSize(iter_is, &num_elem); CHKERRQ(ierr); 214 ierr = ISGetIndices(iter_is, &iter_indices); CHKERRQ(ierr); 215 } else { 216 num_elem = 0; 217 iter_indices = NULL; 218 } 219 ierr = PetscMalloc1(num_elem*PetscPowInt(P, topo_dim), &elem_restr_offsets); 220 CHKERRQ(ierr); 221 for (p = 0, e_offset = 0; p < num_elem; p++) { 222 PetscInt c = iter_indices[p]; 223 PetscInt num_indices, *indices, num_nodes; 224 ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, 225 &num_indices, &indices, NULL, NULL); 226 CHKERRQ(ierr); 227 bool flip = false; 228 if (height > 0) { 229 PetscInt num_cells, num_faces, start = -1; 230 const PetscInt *orients, *faces, *cells; 231 ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr); 232 ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr); 233 if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 234 "Expected one cell in support of exterior face, but got %D cells", 235 num_cells); 236 ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr); 237 ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr); 238 for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;} 239 if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, 240 "Could not find face %D in cone of its support", 241 c); 242 ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr); 243 if (orients[start] < 0) flip = true; 244 } 245 if (num_indices % field_off[num_fields]) SETERRQ1(PETSC_COMM_SELF, 246 PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 247 c); 248 num_nodes = num_indices / field_off[num_fields]; 249 for (PetscInt i = 0; i < num_nodes; i++) { 250 PetscInt ii = i; 251 if (flip) { 252 if (P == num_nodes) ii = num_nodes - 1 - i; 253 else if (P*P == num_nodes) { 254 PetscInt row = i / P, col = i % P; 255 ii = row + col * P; 256 } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, 257 "No support for flipping point with %D nodes != P (%D) or P^2", 258 num_nodes, P); 259 } 260 // Check that indices are blocked by node and thus can be coalesced as a single field with 261 // field_off[num_fields] = sum(num_comp) components. 262 for (PetscInt f = 0; f < num_fields; f++) { 263 for (PetscInt j = 0; j < num_comp[f]; j++) { 264 if (Involute(indices[field_off[f]*num_nodes + ii*num_comp[f] + j]) 265 != Involute(indices[ii*num_comp[0]]) + field_off[f] + j) 266 SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 267 "Cell %D closure indices not interlaced for node %D field %D component %D", 268 c, ii, f, j); 269 } 270 } 271 // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 272 PetscInt loc = Involute(indices[ii*num_comp[0]]); 273 elem_restr_offsets[e_offset++] = loc; 274 } 275 ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, 276 &num_indices, &indices, NULL, NULL); 277 CHKERRQ(ierr); 278 } 279 if (e_offset != num_elem*PetscPowInt(P, topo_dim)) 280 SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB, 281 "ElemRestriction of size (%D,%D) initialized %D nodes", num_elem, 282 PetscPowInt(P, topo_dim),e_offset); 283 if (iter_is) { 284 ierr = ISRestoreIndices(iter_is, &iter_indices); CHKERRQ(ierr); 285 } 286 ierr = ISDestroy(&iter_is); CHKERRQ(ierr); 287 288 ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr); 289 ierr = VecGetLocalSize(U_loc, &num_dof); CHKERRQ(ierr); 290 ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr); 291 CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, topo_dim), 292 field_off[num_fields], 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 293 elem_restr_offsets, elem_restr); 294 ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 }; 297 298 // ----------------------------------------------------------------------------- 299