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 p_degree, PetscInt q_extra, 98 PetscInt num_comp_u, 99 PetscInt dim, bool enforce_bc, BCFunction bc_func) { 100 PetscInt ierr, marker_ids[1] = {1}; 101 PetscInt q_degree = p_degree + q_extra; 102 PetscFE fe; 103 MPI_Comm comm; 104 PetscBool is_simplex = PETSC_TRUE; 105 106 PetscFunctionBeginUser; 107 108 // Check if simplex or tensor-product mesh 109 ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr); 110 // Setup FE 111 ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr); 112 ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree, 113 q_degree, &fe); CHKERRQ(ierr); 114 ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 115 ierr = DMCreateDS(dm); CHKERRQ(ierr); 116 117 { 118 // create FE field for coordinates 119 PetscFE fe_coords; 120 PetscInt num_comp_coord; 121 ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 122 ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree, 123 &fe_coords); CHKERRQ(ierr); 124 ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 125 ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 126 } 127 128 // Setup DM 129 //ierr = DMCreateDS(dm); CHKERRQ(ierr); 130 if (enforce_bc) { 131 PetscBool has_label; 132 DMHasLabel(dm, "marker", &has_label); 133 if (!has_label) {CreateBCLabel(dm, "marker");} 134 DMLabel label; 135 ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr); 136 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, 137 marker_ids, 0, 0, NULL, (void(*)(void))bc_func, 138 NULL, NULL, NULL); CHKERRQ(ierr); 139 } 140 141 if (!is_simplex) { 142 DM dm_coord; 143 ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 144 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 145 CHKERRQ(ierr); 146 ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); 147 CHKERRQ(ierr); 148 } 149 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 150 151 PetscFunctionReturn(0); 152 }; 153 154 // ----------------------------------------------------------------------------- 155 // Get CEED restriction data from DMPlex 156 // ----------------------------------------------------------------------------- 157 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, 158 DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) { 159 PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets; 160 PetscErrorCode ierr; 161 162 PetscFunctionBeginUser; 163 164 ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, 165 &elem_size, &num_comp, &num_dof, &elem_restr_offsets); 166 CHKERRQ(ierr); 167 168 CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 169 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 170 elem_restr_offsets, elem_restr); 171 ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr); 172 173 PetscFunctionReturn(0); 174 }; 175 176 // ----------------------------------------------------------------------------- 177 // Utility function - convert from DMPolytopeType to CeedElemTopology 178 // ----------------------------------------------------------------------------- 179 CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) { 180 switch (cell_type) { 181 case DM_POLYTOPE_TRIANGLE: return CEED_TOPOLOGY_TRIANGLE; 182 case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD; 183 case DM_POLYTOPE_TETRAHEDRON: return CEED_TOPOLOGY_TET; 184 case DM_POLYTOPE_HEXAHEDRON: return CEED_TOPOLOGY_HEX; 185 default: return 0; 186 } 187 } 188 189 // ----------------------------------------------------------------------------- 190 // Get CEED Basis from DMPlex 191 // ----------------------------------------------------------------------------- 192 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, 193 CeedInt label_value, CeedInt height, 194 CeedInt dm_field, CeedBasis *basis) { 195 PetscErrorCode ierr; 196 PetscDS ds; 197 PetscFE fe; 198 PetscQuadrature quadrature; 199 PetscBool is_simplex = PETSC_TRUE; 200 PetscInt dim, ds_field = -1, num_comp, P, Q; 201 202 PetscFunctionBeginUser; 203 204 // Get basis information 205 { 206 IS field_is; 207 const PetscInt *fields; 208 PetscInt num_fields; 209 210 ierr = DMGetRegionDS(dm, domain_label, &field_is, &ds); CHKERRQ(ierr); 211 // Translate dm_field to ds_field 212 ierr = ISGetIndices(field_is, &fields); CHKERRQ(ierr); 213 ierr = ISGetSize(field_is, &num_fields); CHKERRQ(ierr); 214 for (PetscInt i = 0; i < num_fields; i++) { 215 if (dm_field == fields[i]) { 216 ds_field = i; 217 break; 218 } 219 } 220 ierr = ISRestoreIndices(field_is, &fields); CHKERRQ(ierr); 221 } 222 if (ds_field == -1) { 223 // LCOV_EXCL_START 224 SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, 225 "Could not find dm_field %D in DS", dm_field); 226 // LCOV_EXCL_STOP 227 } 228 229 // Get element information 230 { 231 PetscDualSpace dual_space; 232 PetscInt num_dual_basis_vectors; 233 234 ierr = PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe); 235 CHKERRQ(ierr); 236 ierr = PetscFEGetHeightSubspace(fe, height, &fe); CHKERRQ(ierr); 237 ierr = PetscFEGetSpatialDimension(fe, &dim); CHKERRQ(ierr); 238 ierr = PetscFEGetNumComponents(fe, &num_comp); CHKERRQ(ierr); 239 ierr = PetscFEGetDualSpace(fe, &dual_space); CHKERRQ(ierr); 240 ierr = PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors); 241 CHKERRQ(ierr); 242 P = num_dual_basis_vectors / num_comp; 243 ierr = PetscFEGetQuadrature(fe, &quadrature); CHKERRQ(ierr); 244 ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL); 245 CHKERRQ(ierr); 246 } 247 248 // Check if simplex or tensor-product mesh 249 ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr); 250 // Build libCEED basis 251 if (is_simplex) { 252 PetscInt num_derivatives = 1, first_point; 253 PetscInt ids[1] = {label_value}; 254 PetscTabulation basis_tabulation; 255 const PetscScalar *q_points, *q_weights; 256 DMLabel depth_label; 257 DMPolytopeType cell_type; 258 CeedElemTopology elem_topo; 259 PetscScalar *interp, *grad; 260 261 // Use depth label if no domain label present 262 if (!domain_label) { 263 PetscInt depth; 264 265 ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 266 ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr); 267 ids[0] = depth - height; 268 } 269 // Get cell interp, grad, and quadrature data 270 ierr = PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation); 271 CHKERRQ(ierr); 272 ierr = PetscQuadratureGetData(quadrature, NULL, NULL, NULL, &q_points, 273 &q_weights); CHKERRQ(ierr); 274 ierr = DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 275 1, ids, height, &first_point, NULL); 276 CHKERRQ(ierr); 277 ierr = DMPlexGetCellType(dm, first_point, &cell_type); CHKERRQ(ierr); 278 elem_topo = ElemTopologyP2C(cell_type); 279 if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, 280 "DMPlex topology not supported"); 281 // Convert to libCEED orientation 282 ierr = PetscCalloc(P * Q * sizeof(PetscScalar), &interp); CHKERRQ(ierr); 283 ierr = PetscCalloc(P * Q * dim * sizeof(PetscScalar), &grad); CHKERRQ(ierr); 284 const CeedInt c = 0; 285 for (CeedInt q = 0; q < Q; q++) { 286 for (CeedInt p = 0; p < P; p++) { 287 interp[q*P + p] = basis_tabulation->T[0][(q*P + p)*num_comp*num_comp + c]; 288 for (CeedInt d = 0; d < dim; d++) { 289 grad[(d*Q + q)*P + p] = basis_tabulation->T[1][((q*P + p)*num_comp*num_comp + c) 290 *dim + d]; 291 } 292 } 293 } 294 // Finaly, create libCEED basis 295 ierr = CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, 296 q_points, q_weights, basis); 297 CHKERRQ(ierr); 298 ierr = PetscFree(interp); CHKERRQ(ierr); 299 ierr = PetscFree(grad); CHKERRQ(ierr); 300 } else { 301 CeedInt P_1d = (CeedInt) round(pow(P, 1.0 / dim)); 302 CeedInt Q_1d = (CeedInt) round(pow(Q, 1.0 / dim)); 303 304 ierr = CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, 305 CEED_GAUSS, basis); 306 CHKERRQ(ierr); 307 } 308 309 PetscFunctionReturn(0); 310 }; 311 312 // ----------------------------------------------------------------------------- 313 // Utilities 314 // ----------------------------------------------------------------------------- 315 316 // Utility function, compute three factors of an integer 317 static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 318 for (PetscInt d=0, size_left=size; d<3; d++) { 319 PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d))); 320 while (try * (size_left / try) != size_left) try++; 321 m[reverse ? 2-d : d] = try; 322 size_left /= try; 323 } 324 } 325 326 static int Max3(const PetscInt a[3]) { 327 return PetscMax(a[0], PetscMax(a[1], a[2])); 328 } 329 330 static int Min3(const PetscInt a[3]) { 331 return PetscMin(a[0], PetscMin(a[1], a[2])); 332 } 333 334 // ----------------------------------------------------------------------------- 335 // Create distribute dm 336 // ----------------------------------------------------------------------------- 337 PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) { 338 PetscErrorCode ierr; 339 340 PetscFunctionBeginUser; 341 // Setup DM 342 if (rp->read_mesh) { 343 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE, 344 dm); 345 CHKERRQ(ierr); 346 } else { 347 if (rp->user_l_nodes) { 348 // Find a nicely composite number of elements no less than global nodes 349 PetscMPIInt size; 350 ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr); 351 for (PetscInt g_elem = 352 PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim)); 353 ; 354 g_elem++) { 355 Split3(g_elem, rp->mesh_elem, true); 356 if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break; 357 } 358 } 359 ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex, 360 rp->mesh_elem, 361 NULL, NULL, NULL, PETSC_TRUE, dm); CHKERRQ(ierr); 362 } 363 364 ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 365 ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 366 367 PetscFunctionReturn(0); 368 } 369 370 // ----------------------------------------------------------------------------- 371