1 #include "../include/petscutils.h" 2 3 // ----------------------------------------------------------------------------- 4 // Convert PETSc MemType to libCEED MemType 5 // ----------------------------------------------------------------------------- 6 CeedMemType MemTypeP2C(PetscMemType mem_type) { 7 return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; 8 } 9 10 // ----------------------------------------------------------------------------- 11 // Apply 3D Kershaw mesh transformation 12 // ----------------------------------------------------------------------------- 13 // Transition from a value of "a" for x=0, to a value of "b" for x=1. Optionally 14 // smooth -- see the commented versions at the end. 15 static double step(const double a, const double b, double x) { 16 if (x <= 0) return a; 17 if (x >= 1) return b; 18 return a + (b-a) * (x); 19 } 20 21 // 1D transformation at the right boundary 22 static double right(const double eps, const double x) { 23 return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1); 24 } 25 26 // 1D transformation at the left boundary 27 static double left(const double eps, const double x) { 28 return 1-right(eps,1-x); 29 } 30 31 // Apply 3D Kershaw mesh transformation 32 // The eps parameters are in (0, 1] 33 // Uniform mesh is recovered for eps=1 34 PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) { 35 PetscErrorCode ierr; 36 Vec coord; 37 PetscInt ncoord; 38 PetscScalar *c; 39 40 PetscFunctionBeginUser; 41 ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr); 42 ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr); 43 ierr = VecGetArray(coord, &c); CHKERRQ(ierr); 44 45 for (PetscInt i = 0; i < ncoord; i += 3) { 46 PetscScalar x = c[i], y = c[i+1], z = c[i+2]; 47 PetscInt layer = x*6; 48 PetscScalar lambda = (x-layer/6.0)*6; 49 c[i] = x; 50 51 switch (layer) { 52 case 0: 53 c[i+1] = left(eps, y); 54 c[i+2] = left(eps, z); 55 break; 56 case 1: 57 case 4: 58 c[i+1] = step(left(eps, y), right(eps, y), lambda); 59 c[i+2] = step(left(eps, z), right(eps, z), lambda); 60 break; 61 case 2: 62 c[i+1] = step(right(eps, y), left(eps, y), lambda/2); 63 c[i+2] = step(right(eps, z), left(eps, z), lambda/2); 64 break; 65 case 3: 66 c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2); 67 c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2); 68 break; 69 default: 70 c[i+1] = right(eps, y); 71 c[i+2] = right(eps, z); 72 } 73 } 74 ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr); 75 PetscFunctionReturn(0); 76 } 77 78 // ----------------------------------------------------------------------------- 79 // Create BC label 80 // ----------------------------------------------------------------------------- 81 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 82 int ierr; 83 DMLabel label; 84 85 PetscFunctionBeginUser; 86 87 ierr = DMCreateLabel(dm, name); CHKERRQ(ierr); 88 ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr); 89 ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr); 90 91 PetscFunctionReturn(0); 92 }; 93 94 // ----------------------------------------------------------------------------- 95 // This function sets up a DM for a given degree 96 // ----------------------------------------------------------------------------- 97 PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u, 98 PetscInt dim, bool enforce_bc, BCFunction bc_func) { 99 PetscInt ierr, marker_ids[1] = {1}; 100 PetscFE fe; 101 MPI_Comm comm; 102 PetscBool is_simplex = PETSC_TRUE; 103 104 PetscFunctionBeginUser; 105 106 // Check if simplex or tensor-product mesh 107 ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr); 108 // Setup FE 109 ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr); 110 ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, degree, degree, 111 &fe); CHKERRQ(ierr); 112 ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 113 ierr = DMCreateDS(dm); CHKERRQ(ierr); 114 { 115 DM dm_coord; 116 PetscDS ds_coord; 117 PetscFE fe_coord_current, fe_coord_new; 118 PetscDualSpace fe_coord_dual_space; 119 PetscInt fe_coord_order, num_comp_coord; 120 121 ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 122 ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 123 ierr = DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord); CHKERRQ(ierr); 124 ierr = PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current); CHKERRQ(ierr); 125 ierr = PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space); CHKERRQ(ierr); 126 ierr = PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order); CHKERRQ(ierr); 127 128 // Create FE for coordinates 129 ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, degree, &fe_coord_new); 130 CHKERRQ(ierr); 131 ierr = DMProjectCoordinates(dm, fe_coord_new); CHKERRQ(ierr); 132 ierr = PetscFEDestroy(&fe_coord_new); CHKERRQ(ierr); 133 } 134 /* 135 { 136 // create FE field for coordinates 137 PetscFE fe_coords; 138 PetscInt num_comp_coord; 139 ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 140 ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, degree, 141 &fe_coords); CHKERRQ(ierr); 142 ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 143 ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 144 } 145 */ 146 // Setup DM 147 //ierr = DMCreateDS(dm); CHKERRQ(ierr); 148 if (enforce_bc) { 149 PetscBool has_label; 150 DMHasLabel(dm, "marker", &has_label); 151 if (!has_label) {CreateBCLabel(dm, "marker");} 152 DMLabel label; 153 ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr); 154 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, 155 marker_ids, 0, 0, NULL, (void(*)(void))bc_func, 156 NULL, NULL, NULL); CHKERRQ(ierr); 157 } 158 159 if (!is_simplex) { 160 DM dm_coord; 161 ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 162 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); CHKERRQ(ierr); 163 ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); CHKERRQ(ierr); 164 } 165 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 166 167 PetscFunctionReturn(0); 168 }; 169 170 // ----------------------------------------------------------------------------- 171 // Utility function - essential BC dofs are encoded in closure indices as -(i+1) 172 // ----------------------------------------------------------------------------- 173 PetscInt Involute(PetscInt i) { 174 return i >= 0 ? i : -(i + 1); 175 }; 176 177 // ----------------------------------------------------------------------------- 178 // Get CEED restriction data from DMPlex 179 // ----------------------------------------------------------------------------- 180 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, 181 DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) { 182 PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets; 183 PetscErrorCode ierr; 184 185 PetscFunctionBeginUser; 186 187 ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, 188 &elem_size, &num_comp, &num_dof, &elem_restr_offsets); 189 CHKERRQ(ierr); 190 191 CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 192 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 193 elem_restr_offsets, elem_restr); 194 ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr); 195 196 PetscFunctionReturn(0); 197 }; 198 199 // ----------------------------------------------------------------------------- 200 // Utility function - convert from DMPolytopeType to CeedElemTopology 201 // ----------------------------------------------------------------------------- 202 static inline CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) { 203 switch (cell_type) { 204 case DM_POLYTOPE_TRIANGLE: return CEED_TOPOLOGY_TRIANGLE; 205 case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD; 206 case DM_POLYTOPE_TETRAHEDRON: return CEED_TOPOLOGY_TET; 207 case DM_POLYTOPE_HEXAHEDRON: return CEED_TOPOLOGY_HEX; 208 default: return 0; 209 } 210 } 211 212 // ----------------------------------------------------------------------------- 213 // Get CEED Basis from DMPlex 214 // ----------------------------------------------------------------------------- 215 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, 216 CeedInt label_value, CeedInt height, 217 CeedInt dm_field, CeedBasis *basis) { 218 PetscErrorCode ierr; 219 PetscDS ds; 220 PetscFE fe; 221 PetscQuadrature quadrature; 222 PetscBool is_simplex = PETSC_TRUE; 223 PetscInt dim, ds_field = -1, num_comp, P, Q; 224 225 PetscFunctionBeginUser; 226 227 // Get basis information 228 { 229 IS field_is; 230 const PetscInt *fields; 231 PetscInt num_fields; 232 233 ierr = DMGetRegionDS(dm, domain_label, &field_is, &ds); CHKERRQ(ierr); 234 // Translate dm_field to ds_field 235 ierr = ISGetIndices(field_is, &fields); CHKERRQ(ierr); 236 ierr = ISGetSize(field_is, &num_fields); CHKERRQ(ierr); 237 for (PetscInt i = 0; i < num_fields; i++) { 238 if (dm_field == fields[i]) { 239 ds_field = i; 240 break; 241 } 242 } 243 ierr = ISRestoreIndices(field_is, &fields); CHKERRQ(ierr); 244 } 245 if (ds_field == -1) { 246 // LCOV_EXCL_START 247 SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, 248 "Could not find dm_field %D in DS", dm_field); 249 // LCOV_EXCL_STOP 250 } 251 252 // Get element information 253 { 254 PetscDualSpace dual_space; 255 PetscInt num_dual_basis_vectors; 256 257 ierr = PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe); 258 CHKERRQ(ierr); 259 ierr = PetscFEGetHeightSubspace(fe, height, &fe); CHKERRQ(ierr); 260 ierr = PetscFEGetSpatialDimension(fe, &dim); CHKERRQ(ierr); 261 ierr = PetscFEGetNumComponents(fe, &num_comp); CHKERRQ(ierr); 262 ierr = PetscFEGetDualSpace(fe, &dual_space); CHKERRQ(ierr); 263 ierr = PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors); 264 CHKERRQ(ierr); 265 P = num_dual_basis_vectors / num_comp; 266 ierr = PetscFEGetQuadrature(fe, &quadrature); CHKERRQ(ierr); 267 ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL); 268 CHKERRQ(ierr); 269 } 270 271 // Check if simplex or tensor-product mesh 272 ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr); 273 // Build libCEED basis 274 if (is_simplex) { 275 PetscInt num_derivatives = 1, first_point; 276 PetscInt ids[1] = {label_value}; 277 PetscTabulation basis_tabulation; 278 const PetscScalar *q_points, *q_weights; 279 DMLabel depth_label; 280 DMPolytopeType cell_type; 281 CeedElemTopology elem_topo; 282 PetscScalar *interp, *grad; 283 284 // Use depth label if no domain label present 285 if (!domain_label) { 286 PetscInt depth; 287 288 ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 289 ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr); 290 ids[0] = depth - height; 291 } 292 // Get cell interp, grad, and quadrature data 293 ierr = PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation); 294 CHKERRQ(ierr); 295 ierr = PetscQuadratureGetData(quadrature, NULL, NULL, NULL, &q_points, 296 &q_weights); CHKERRQ(ierr); 297 ierr = DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 298 1, ids, height, &first_point, NULL); 299 CHKERRQ(ierr); 300 ierr = DMPlexGetCellType(dm, first_point, &cell_type); CHKERRQ(ierr); 301 elem_topo = ElemTopologyP2C(cell_type); 302 if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, 303 "DMPlex topology not supported"); 304 // Convert to libCEED orientation 305 ierr = PetscCalloc(P * Q * sizeof(PetscScalar), &interp); CHKERRQ(ierr); 306 ierr = PetscCalloc(P * Q * dim * sizeof(PetscScalar), &grad); CHKERRQ(ierr); 307 const CeedInt c = 0; 308 for (CeedInt q = 0; q < Q; q++) { 309 for (CeedInt p = 0; p < P; p++) { 310 interp[q*P + p] = basis_tabulation->T[0][(q*P + p)*num_comp*num_comp + c]; 311 for (CeedInt d = 0; d < dim; d++) { 312 grad[(d*Q + q)*P + p] = basis_tabulation->T[1][((q*P + p)*num_comp*num_comp + c) 313 *dim + d]; 314 } 315 } 316 } 317 // Finaly, create libCEED basis 318 ierr = CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, 319 q_points, q_weights, basis); 320 CHKERRQ(ierr); 321 ierr = PetscFree(interp); CHKERRQ(ierr); 322 ierr = PetscFree(grad); CHKERRQ(ierr); 323 } else { 324 CeedInt P_1d = (CeedInt) round(pow(P, 1.0 / dim)); 325 CeedInt Q_1d = (CeedInt) round(pow(Q, 1.0 / dim)); 326 327 ierr = CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, 328 CEED_GAUSS, basis); 329 CHKERRQ(ierr); 330 } 331 332 PetscFunctionReturn(0); 333 }; 334 335 // ----------------------------------------------------------------------------- 336