1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utilities for setting up DM and PetscFE 6 7 #include <ceed.h> 8 #include <dm-utils.h> 9 #include <petsc-ceed-utils.h> 10 #include <petscdmplex.h> 11 #include <petscds.h> 12 13 /** 14 @brief Convert `DM` field to `DS` field. 15 16 @param[in] dm `DM` holding mesh 17 @param[in] domain_label Label for `DM` domain 18 @param[in] dm_field Index of `DM` field 19 @param[out] ds_field Index of `DS` field 20 21 @return An error code: 0 - success, otherwise - failure 22 **/ 23 PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) { 24 PetscDS ds; 25 IS field_is; 26 const PetscInt *fields; 27 PetscInt num_fields; 28 29 PetscFunctionBeginUser; 30 PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL)); 31 PetscCall(ISGetIndices(field_is, &fields)); 32 PetscCall(ISGetSize(field_is, &num_fields)); 33 for (PetscInt i = 0; i < num_fields; i++) { 34 if (dm_field == fields[i]) { 35 *ds_field = i; 36 break; 37 } 38 } 39 PetscCall(ISRestoreIndices(field_is, &fields)); 40 41 PetscCheck(*ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 42 PetscFunctionReturn(PETSC_SUCCESS); 43 } 44 45 /** 46 @brief Destroy `ElemRestriction` in a `PetscContainer`. 47 48 Not collective across MPI processes. 49 50 @param[in,out] ctx `CeedElemRestriction` 51 52 @return An error code: 0 - success, otherwise - failure 53 **/ 54 static PetscErrorCode DMPlexCeedElemRestrictionDestroy(void **ctx) { 55 CeedElemRestriction rstr = *(CeedElemRestriction *)ctx; 56 57 PetscFunctionBegin; 58 PetscCallCeed(CeedElemRestrictionReturnCeed(rstr), CeedElemRestrictionDestroy(&rstr)); 59 PetscFunctionReturn(PETSC_SUCCESS); 60 } 61 62 /** 63 @brief Create `CeedElemRestriction` from `DMPlex`. 64 65 The `CeedElemRestriction` is stored in the `DMPlex` object for reuse; 66 if the routine is called again with the same arguments, the same `CeedElemRestriction` is returned. 67 68 Not collective across MPI processes. 69 70 @param[in] ceed `Ceed` context 71 @param[in] dm `DMPlex` holding mesh 72 @param[in] domain_label `DMLabel` for `DMPlex` domain 73 @param[in] label_value Stratum value 74 @param[in] height Height of `DMPlex` topology 75 @param[in] dm_field Index of `DMPlex` field 76 @param[out] restriction `CeedElemRestriction` for `DMPlex` 77 78 @return An error code: 0 - success, otherwise - failure 79 **/ 80 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 81 CeedElemRestriction *restriction) { 82 char container_name[1024]; 83 CeedElemRestriction container_restriction; 84 85 PetscFunctionBeginUser; 86 { 87 const char *label_name = NULL; 88 89 if (domain_label) PetscCall(PetscObjectGetName((PetscObject)domain_label, &label_name)); 90 PetscCall(PetscSNPrintf(container_name, sizeof(container_name), 91 "DM CeedElemRestriction - label: %s label value: %" PetscInt_FMT " height: %" PetscInt_FMT " DM field: %" PetscInt_FMT, 92 label_name ? label_name : "NONE", label_value, height, dm_field)); 93 } 94 PetscCall(PetscObjectContainerQuery((PetscObject)dm, container_name, (void **)&container_restriction)); 95 96 if (container_restriction) { 97 // Copy existing restriction 98 *restriction = NULL; 99 PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(container_restriction, restriction)); 100 } else { 101 PetscInt num_elem, elem_size, num_dof, num_comp, *restriction_offsets_petsc; 102 CeedInt *restriction_offsets_ceed = NULL; 103 104 // Build restriction 105 PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, 106 &restriction_offsets_petsc)); 107 PetscCall(IntArrayPetscToCeed(num_elem * elem_size, &restriction_offsets_petsc, &restriction_offsets_ceed)); 108 PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 109 restriction_offsets_ceed, restriction)); 110 PetscCall(PetscFree(restriction_offsets_ceed)); 111 112 // Set in container 113 PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(*restriction, &container_restriction)); 114 PetscCall(PetscObjectContainerCompose((PetscObject)dm, container_name, container_restriction, DMPlexCeedElemRestrictionDestroy)); 115 } 116 PetscFunctionReturn(PETSC_SUCCESS); 117 } 118 119 /** 120 @brief Create `CeedElemRestriction` from `DMPlex` domain for mesh coordinates. 121 122 The `CeedElemRestriction` is stored in the coordinate `DMPlex` object for reuse; 123 if the routine is called again with the same arguments, the same `CeedElemRestriction` is returned. 124 125 Not collective across MPI processes. 126 127 @param[in] ceed `Ceed` context 128 @param[in] dm `DMPlex` holding mesh 129 @param[in] domain_label Label for `DMPlex` domain 130 @param[in] label_value Stratum value 131 @param[in] height Height of `DMPlex` topology 132 @param[out] restriction `CeedElemRestriction` for mesh 133 134 @return An error code: 0 - success, otherwise - failure 135 **/ 136 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 137 CeedElemRestriction *restriction) { 138 DM dm_coord; 139 140 PetscFunctionBeginUser; 141 PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 142 if (!dm_coord) { 143 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 144 } 145 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm_coord, domain_label, label_value, height, 0, restriction)); 146 PetscFunctionReturn(PETSC_SUCCESS); 147 } 148 149 /** 150 @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data. 151 152 Not collective across MPI processes. 153 154 @param[in] ceed `Ceed` context 155 @param[in] dm `DMPlex` holding mesh 156 @param[in] domain_label Label for `DMPlex` domain 157 @param[in] label_value Stratum value 158 @param[in] height Height of `DMPlex` topology 159 @param[in] q_data_size Number of components for `QFunction` data 160 @param[in] is_collocated Boolean flag indicating if the data is collocated on the nodes (`PETSC_TRUE`) or on quadrature points (`PETSC_FALSE`) 161 @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 162 163 @return An error code: 0 - success, otherwise - failure 164 **/ 165 static PetscErrorCode DMPlexCeedElemRestrictionStridedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 166 PetscInt q_data_size, PetscBool is_collocated, CeedElemRestriction *restriction) { 167 PetscInt num_elem, num_qpts, dm_field = 0; 168 169 PetscFunctionBeginUser; 170 { // Get number of elements 171 PetscInt depth; 172 DMLabel depth_label; 173 IS point_is, depth_is; 174 175 PetscCall(DMPlexGetDepth(dm, &depth)); 176 PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 177 PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 178 if (domain_label) { 179 IS domain_is; 180 181 PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is)); 182 if (domain_is) { 183 PetscCall(ISIntersect(depth_is, domain_is, &point_is)); 184 PetscCall(ISDestroy(&domain_is)); 185 } else { 186 point_is = NULL; 187 } 188 PetscCall(ISDestroy(&depth_is)); 189 } else { 190 point_is = depth_is; 191 } 192 if (point_is) { 193 PetscCall(ISGetLocalSize(point_is, &num_elem)); 194 } else { 195 num_elem = 0; 196 } 197 PetscCall(ISDestroy(&point_is)); 198 } 199 200 { // Get number of quadrature points 201 PetscDS ds; 202 PetscFE fe; 203 PetscInt ds_field = -1; 204 205 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 206 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 207 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 208 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 209 if (is_collocated) { 210 PetscDualSpace dual_space; 211 PetscInt num_dual_basis_vectors, dim, num_comp; 212 213 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 214 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 215 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 216 PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 217 num_qpts = num_dual_basis_vectors / num_comp; 218 } else { 219 PetscQuadrature quadrature; 220 221 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 222 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 223 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 224 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 225 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 226 PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &num_qpts, NULL, NULL)); 227 } 228 } 229 230 // Create the restriction 231 PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, num_elem * num_qpts * q_data_size, CEED_STRIDES_BACKEND, 232 restriction)); 233 PetscFunctionReturn(PETSC_SUCCESS); 234 } 235 236 /** 237 @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data. 238 239 Not collective across MPI processes. 240 241 @param[in] ceed `Ceed` context 242 @param[in] dm `DMPlex` holding mesh 243 @param[in] domain_label Label for `DMPlex` domain 244 @param[in] label_value Stratum value 245 @param[in] height Height of `DMPlex` topology 246 @param[in] q_data_size Number of components for `QFunction` data 247 @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 248 249 @return An error code: 0 - success, otherwise - failure 250 **/ 251 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 252 PetscInt q_data_size, CeedElemRestriction *restriction) { 253 PetscFunctionBeginUser; 254 PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_FALSE, restriction)); 255 PetscFunctionReturn(PETSC_SUCCESS); 256 } 257 258 /** 259 @brief Create `CeedElemRestriction` from `DMPlex` domain for nodally collocated auxilury `QFunction` data. 260 261 Not collective across MPI processes. 262 263 @param[in] ceed `Ceed` context 264 @param[in] dm `DMPlex` holding mesh 265 @param[in] domain_label Label for `DMPlex` domain 266 @param[in] label_value Stratum value 267 @param[in] height Height of `DMPlex` topology 268 @param[in] q_data_size Number of components for `QFunction` data 269 @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 270 271 @return An error code: 0 - success, otherwise - failure 272 **/ 273 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 274 PetscInt q_data_size, CeedElemRestriction *restriction) { 275 PetscFunctionBeginUser; 276 PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_TRUE, restriction)); 277 PetscFunctionReturn(PETSC_SUCCESS); 278 } 279 280 /** 281 @brief Creates a tensor-product uniform quadrature. 282 This is only for comparison studies and generally should not be used. 283 284 @param[in] dim Spatial dimension 285 @param[in] num_comp Number of components 286 @param[in] num_points Number of points in one dimension 287 @param[in] a Left end of interval (often -1) 288 @param[in] b Right end of interval (often +1) 289 @param[out] q `PetscQuadrature` object 290 291 @return An error code: 0 - success, otherwise - failure 292 **/ 293 static PetscErrorCode PetscDTUniformTensorQuadrature(PetscInt dim, PetscInt num_comp, PetscInt num_points, PetscReal a, PetscReal b, 294 PetscQuadrature *q) { 295 DMPolytopeType ct; 296 PetscInt num_total_points = dim > 1 ? (dim > 2 ? (num_points * num_points * num_points) : (num_points * num_points)) : num_points; 297 PetscReal *coords, *weights, *coords_1d; 298 299 PetscFunctionBeginUser; 300 PetscCall(PetscMalloc1(num_total_points * dim, &coords)); 301 PetscCall(PetscMalloc1(num_total_points * num_comp, &weights)); 302 // Compute weights, points 303 switch (dim) { 304 case 0: 305 ct = DM_POLYTOPE_POINT; 306 PetscCall(PetscFree(coords)); 307 PetscCall(PetscFree(weights)); 308 PetscCall(PetscMalloc1(1, &coords)); 309 PetscCall(PetscMalloc1(num_comp, &weights)); 310 num_total_points = 1; 311 coords[0] = 0.0; 312 for (PetscInt c = 0; c < num_comp; c++) weights[c] = 1.0; 313 break; 314 case 1: { 315 ct = DM_POLYTOPE_SEGMENT; 316 PetscReal step = (b - a) / num_points; 317 318 for (PetscInt i = 0; i < num_points; i++) { 319 coords[i] = step * (i + 0.5) + a; 320 for (PetscInt c = 0; c < num_comp; c++) weights[i * num_comp + c] = 1.0; 321 } 322 } break; 323 case 2: { 324 ct = DM_POLYTOPE_QUADRILATERAL; 325 PetscCall(PetscMalloc1(num_points, &coords_1d)); 326 PetscReal step = (b - a) / num_points; 327 328 for (PetscInt i = 0; i < num_points; i++) coords_1d[i] = step * (i + 0.5) + a; 329 for (PetscInt i = 0; i < num_points; i++) { 330 for (PetscInt j = 0; j < num_points; j++) { 331 coords[(i * num_points + j) * dim + 0] = coords_1d[i]; 332 coords[(i * num_points + j) * dim + 1] = coords_1d[j]; 333 for (PetscInt c = 0; c < num_comp; c++) weights[(i * num_points + j) * num_comp + c] = 1.0; 334 } 335 } 336 PetscCall(PetscFree(coords_1d)); 337 } break; 338 case 3: { 339 ct = DM_POLYTOPE_HEXAHEDRON; 340 PetscCall(PetscMalloc1(num_points, &coords_1d)); 341 PetscReal step = (b - a) / num_points; 342 343 for (PetscInt i = 0; i < num_points; i++) coords_1d[i] = step * (i + 0.5) + a; 344 for (PetscInt i = 0; i < num_points; i++) { 345 for (PetscInt j = 0; j < num_points; j++) { 346 for (PetscInt k = 0; k < num_points; k++) { 347 coords[((i * num_points + j) * num_points + k) * dim + 0] = coords_1d[i]; 348 coords[((i * num_points + j) * num_points + k) * dim + 1] = coords_1d[j]; 349 coords[((i * num_points + j) * num_points + k) * dim + 2] = coords_1d[k]; 350 for (PetscInt c = 0; c < num_comp; c++) weights[((i * num_points + j) * num_points + k) * num_comp + c] = 1.0; 351 } 352 } 353 } 354 PetscCall(PetscFree(coords_1d)); 355 } break; 356 // LCOV_EXCL_START 357 default: 358 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %" PetscInt_FMT, dim); 359 // LCOV_EXCL_STOP 360 } 361 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q)); 362 PetscCall(PetscQuadratureSetCellType(*q, ct)); 363 PetscCall(PetscQuadratureSetOrder(*q, num_points)); 364 PetscCall(PetscQuadratureSetData(*q, dim, num_comp, num_total_points, coords, weights)); 365 PetscCall(PetscObjectChangeTypeName((PetscObject)*q, "UniformTensor")); 366 PetscFunctionReturn(PETSC_SUCCESS); 367 } 368 369 /** 370 @brief Setup `DM` with FE space of appropriate degree 371 372 @param[in] comm MPI communicator 373 @param[in] dim Spatial dimension 374 @param[in] num_comp Number of components 375 @param[in] is_simplex Flag for simplex or tensor product element 376 @param[in] order Polynomial order of space 377 @param[in] q_order Quadrature order 378 @param[in] prefix The options prefix, or `NULL` 379 @param[out] fem `PetscFE` object 380 381 @return An error code: 0 - success, otherwise - failure 382 **/ 383 PetscErrorCode PetscFECreateLagrangeFromOptions(MPI_Comm comm, PetscInt dim, PetscInt num_comp, PetscBool is_simplex, PetscInt order, 384 PetscInt q_order, const char prefix[], PetscFE *fem) { 385 PetscBool is_tensor = !is_simplex; 386 DMPolytopeType polytope_type = DMPolytopeTypeSimpleShape(dim, is_simplex); 387 PetscSpace fe_space; 388 PetscDualSpace fe_dual_space; 389 PetscQuadrature quadrature, face_quadrature; 390 391 PetscFunctionBeginUser; 392 // Create space 393 PetscCall(PetscSpaceCreate(comm, &fe_space)); 394 PetscCall(PetscSpaceSetType(fe_space, PETSCSPACEPOLYNOMIAL)); 395 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fe_space, prefix)); 396 PetscCall(PetscSpacePolynomialSetTensor(fe_space, is_tensor)); 397 PetscCall(PetscSpaceSetNumComponents(fe_space, num_comp)); 398 PetscCall(PetscSpaceSetNumVariables(fe_space, dim)); 399 if (order >= 0) { 400 PetscCall(PetscSpaceSetDegree(fe_space, order, PETSC_DETERMINE)); 401 if (polytope_type == DM_POLYTOPE_TRI_PRISM || polytope_type == DM_POLYTOPE_TRI_PRISM_TENSOR) { 402 PetscSpace fe_space_end, fe_space_side; 403 404 PetscCall(PetscSpaceSetNumComponents(fe_space, 1)); 405 PetscCall(PetscSpaceCreate(comm, &fe_space_end)); 406 PetscCall(PetscSpaceSetType(fe_space_end, PETSCSPACEPOLYNOMIAL)); 407 PetscCall(PetscSpacePolynomialSetTensor(fe_space_end, PETSC_FALSE)); 408 PetscCall(PetscSpaceSetNumComponents(fe_space_end, 1)); 409 PetscCall(PetscSpaceSetNumVariables(fe_space_end, dim - 1)); 410 PetscCall(PetscSpaceSetDegree(fe_space_end, order, PETSC_DETERMINE)); 411 PetscCall(PetscSpaceCreate(comm, &fe_space_side)); 412 PetscCall(PetscSpaceSetType(fe_space_side, PETSCSPACEPOLYNOMIAL)); 413 PetscCall(PetscSpacePolynomialSetTensor(fe_space_side, PETSC_FALSE)); 414 PetscCall(PetscSpaceSetNumComponents(fe_space_side, 1)); 415 PetscCall(PetscSpaceSetNumVariables(fe_space_side, 1)); 416 PetscCall(PetscSpaceSetDegree(fe_space_side, order, PETSC_DETERMINE)); 417 PetscCall(PetscSpaceSetType(fe_space, PETSCSPACETENSOR)); 418 PetscCall(PetscSpaceTensorSetNumSubspaces(fe_space, 2)); 419 PetscCall(PetscSpaceTensorSetSubspace(fe_space, 0, fe_space_end)); 420 PetscCall(PetscSpaceTensorSetSubspace(fe_space, 1, fe_space_side)); 421 PetscCall(PetscSpaceDestroy(&fe_space_end)); 422 PetscCall(PetscSpaceDestroy(&fe_space_side)); 423 424 if (num_comp > 1) { 425 PetscSpace scalar_fe_space = fe_space; 426 427 PetscCall(PetscSpaceCreate(comm, &fe_space)); 428 PetscCall(PetscSpaceSetNumVariables(fe_space, dim)); 429 PetscCall(PetscSpaceSetNumComponents(fe_space, num_comp)); 430 PetscCall(PetscSpaceSetType(fe_space, PETSCSPACESUM)); 431 PetscCall(PetscSpaceSumSetNumSubspaces(fe_space, num_comp)); 432 PetscCall(PetscSpaceSumSetConcatenate(fe_space, PETSC_TRUE)); 433 PetscCall(PetscSpaceSumSetInterleave(fe_space, PETSC_TRUE, PETSC_FALSE)); 434 for (PetscInt i = 0; i < num_comp; i++) PetscCall(PetscSpaceSumSetSubspace(fe_space, i, scalar_fe_space)); 435 PetscCall(PetscSpaceDestroy(&scalar_fe_space)); 436 } 437 } 438 } 439 PetscCall(PetscSpaceSetFromOptions(fe_space)); 440 PetscCall(PetscSpaceSetUp(fe_space)); 441 PetscCall(PetscSpaceGetDegree(fe_space, &order, NULL)); 442 PetscCall(PetscSpacePolynomialGetTensor(fe_space, &is_tensor)); 443 PetscCall(PetscSpaceGetNumComponents(fe_space, &num_comp)); 444 445 // Create dual space 446 PetscCall(PetscDualSpaceCreate(comm, &fe_dual_space)); 447 PetscCall(PetscDualSpaceSetType(fe_dual_space, PETSCDUALSPACELAGRANGE)); 448 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fe_dual_space, prefix)); 449 { 450 DM dual_space_dm; 451 452 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, polytope_type, &dual_space_dm)); 453 PetscCall(PetscDualSpaceSetDM(fe_dual_space, dual_space_dm)); 454 PetscCall(DMDestroy(&dual_space_dm)); 455 } 456 PetscCall(PetscDualSpaceSetNumComponents(fe_dual_space, num_comp)); 457 PetscCall(PetscDualSpaceSetOrder(fe_dual_space, order)); 458 PetscCall(PetscDualSpaceLagrangeSetTensor(fe_dual_space, (is_tensor || (polytope_type == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE)); 459 PetscCall(PetscDualSpaceSetFromOptions(fe_dual_space)); 460 PetscCall(PetscDualSpaceSetUp(fe_dual_space)); 461 462 // Create quadrature 463 q_order = q_order >= 0 ? q_order : order; 464 PetscBool use_uniform = PETSC_FALSE; 465 466 PetscOptionsBegin(comm, NULL, "Uniform quadrature check", NULL); 467 PetscCall(PetscOptionsBool("-petscdt_uniform_tensor_quadrature", "", NULL, use_uniform, &use_uniform, NULL)); 468 PetscOptionsEnd(); 469 PetscCheck(!use_uniform || (is_tensor || dim == 1), comm, PETSC_ERR_SUP, "Can only use uniform quadrature on tensor product elements"); 470 if (use_uniform) { 471 PetscCall(PetscDTUniformTensorQuadrature(dim, 1, q_order + 1, -1.0, 1.0, &quadrature)); 472 PetscCall(PetscDTUniformTensorQuadrature(dim - 1, 1, q_order + 1, -1.0, 1.0, &face_quadrature)); 473 } else { 474 PetscCall(PetscDTCreateDefaultQuadrature(polytope_type, q_order, &quadrature, &face_quadrature)); 475 } 476 477 // Create finite element 478 PetscCall(PetscFECreateFromSpaces(fe_space, fe_dual_space, quadrature, face_quadrature, fem)); 479 PetscCall(PetscFESetFromOptions(*fem)); 480 PetscFunctionReturn(PETSC_SUCCESS); 481 } 482 483 /** 484 @brief Get global `DMPlex` topology type. 485 486 Collective across MPI processes. 487 488 @param[in] dm `DMPlex` holding mesh 489 @param[in] domain_label `DMLabel` for `DMPlex` domain 490 @param[in] label_value Stratum value 491 @param[in] height Height of `DMPlex` topology 492 @param[out] cell_type `DMPlex` topology type 493 **/ 494 static inline PetscErrorCode GetGlobalDMPlexPolytopeType(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 495 DMPolytopeType *cell_type) { 496 PetscInt first_point; 497 PetscInt ids[1] = {label_value}; 498 DMLabel depth_label; 499 500 PetscFunctionBeginUser; 501 // Use depth label if no domain label present 502 if (!domain_label) { 503 PetscInt depth; 504 505 PetscCall(DMPlexGetDepth(dm, &depth)); 506 PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 507 ids[0] = depth - height; 508 } 509 510 // Get cell interp, grad, and quadrature data 511 PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL)); 512 if (first_point != -1) PetscCall(DMPlexGetCellType(dm, first_point, cell_type)); 513 514 { // Form agreement about CellType 515 PetscInt cell_type_local = -1, cell_type_global = -1; 516 517 if (first_point != -1) cell_type_local = (PetscInt)*cell_type; 518 PetscCallMPI(MPIU_Allreduce(&cell_type_local, &cell_type_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 519 *cell_type = (DMPolytopeType)cell_type_global; 520 } 521 PetscFunctionReturn(PETSC_SUCCESS); 522 } 523 524 /** 525 @brief Get the permutation and field offset for a given depth. 526 527 Not collective across MPI processes. 528 529 @param[in] dm `DMPlex` holding mesh 530 @param[in] depth Depth of `DMPlex` topology 531 @param[in] field Index of `DMPlex` field 532 @param[out] permutation Permutation for `DMPlex` field 533 @param[out] field_offset Offset to `DMPlex` field 534 **/ 535 static inline PetscErrorCode GetClosurePermutationAndFieldOffsetAtDepth(DM dm, PetscInt depth, PetscInt field, IS *permutation, 536 PetscInt *field_offset) { 537 PetscBool is_field_continuous = PETSC_TRUE; 538 PetscInt dim, num_fields, offset = 0, size = 0; 539 PetscSection section; 540 541 PetscFunctionBeginUser; 542 *field_offset = 0; 543 544 PetscCall(DMGetDimension(dm, &dim)); 545 PetscCall(DMGetLocalSection(dm, §ion)); 546 PetscCall(PetscSectionGetNumFields(section, &num_fields)); 547 // ---- Permutation size and offset to current field 548 for (PetscInt f = 0; f < num_fields; f++) { 549 PetscBool is_continuous; 550 PetscInt num_components, num_dof_1d, dual_space_size, field_size; 551 PetscObject obj; 552 PetscFE fe = NULL; 553 PetscDualSpace dual_space; 554 555 PetscCall(PetscSectionGetFieldComponents(section, f, &num_components)); 556 PetscCall(DMGetField(dm, f, NULL, &obj)); 557 fe = (PetscFE)obj; 558 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 559 PetscCall(PetscDualSpaceLagrangeGetContinuity(dual_space, &is_continuous)); 560 if (f == field) is_field_continuous = is_continuous; 561 if (!is_continuous) continue; 562 PetscCall(PetscDualSpaceGetDimension(dual_space, &dual_space_size)); 563 num_dof_1d = (PetscInt)PetscCeilReal(PetscPowReal(dual_space_size / num_components, 1.0 / dim)); 564 field_size = PetscPowInt(num_dof_1d, depth) * num_components; 565 if (f < field) offset += field_size; 566 size += field_size; 567 } 568 569 if (is_field_continuous) { 570 *field_offset = offset; 571 PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, depth, size, permutation)); 572 } 573 PetscFunctionReturn(PETSC_SUCCESS); 574 } 575 576 /** 577 @brief Compute field tabulation from `PetscTabulation`. 578 579 Not collective across MPI processes. 580 581 @param[in] dm `DMPlex` holding mesh 582 @param[in] field Index of `DMPlex` field 583 @param[in] face Index of basis face, or 0 584 @param[in] tabulation PETSc basis tabulation 585 @param[out] interp Interpolation matrix in libCEED orientation 586 @param[out] grad Gradient matrix in libCEED orientation 587 588 @note `interp` and `grad` are allocated using `PetscCalloc1` and must be freed by the user. 589 **/ 590 static inline PetscErrorCode ComputeFieldTabulationP2C(DM dm, PetscInt field, PetscInt face, PetscTabulation tabulation, CeedScalar **interp, 591 CeedScalar **grad) { 592 PetscBool is_simplex, has_permutation = PETSC_FALSE; 593 PetscInt field_offset = 0, num_comp, P, Q, dim; 594 const PetscInt *permutation_indices; 595 IS permutation = NULL; 596 597 PetscFunctionBeginUser; 598 num_comp = tabulation->Nc; 599 P = tabulation->Nb / num_comp; 600 Q = tabulation->Np; 601 dim = tabulation->cdim; 602 603 PetscCall(PetscCalloc1(P * Q, interp)); 604 PetscCall(PetscCalloc1(P * Q * dim, grad)); 605 606 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 607 if (!is_simplex) { 608 PetscCall(GetClosurePermutationAndFieldOffsetAtDepth(dm, dim, field, &permutation, &field_offset)); 609 has_permutation = !!permutation; 610 if (has_permutation) PetscCall(ISGetIndices(permutation, &permutation_indices)); 611 } 612 613 const CeedInt c = 0; 614 for (CeedInt q = 0; q < Q; q++) { 615 for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) { 616 CeedInt p_petsc = has_permutation ? (permutation_indices[field_offset + p_ceed * num_comp] - field_offset) : (p_ceed * num_comp); 617 (*interp)[q * P + p_ceed] = tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c]; 618 for (CeedInt d = 0; d < dim; d++) { 619 (*grad)[(d * Q + q) * P + p_ceed] = tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d]; 620 } 621 } 622 } 623 624 if (has_permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices)); 625 PetscCall(ISDestroy(&permutation)); 626 PetscFunctionReturn(PETSC_SUCCESS); 627 } 628 629 /** 630 @brief Get quadrature data from `PetscQuadrature` for use with libCEED. 631 632 Not collective across MPI processes. 633 634 @param[in] fe `PetscFE` object 635 @param[in] quadrature PETSc basis quadrature 636 @param[out] q_dim Quadrature dimension, or NULL if not needed 637 @param[out] num_comp Number of components, or NULL if not needed 638 @param[out] Q Number of quadrature points, or NULL if not needed 639 @param[out] q_points Quadrature points in libCEED orientation, or NULL if not needed 640 @param[out] q_weights Quadrature weights in libCEED orientation, or NULL if not needed 641 642 @note `q_weights` and `q_points` are allocated using `PetscCalloc1` and must be freed by the user. 643 **/ 644 static inline PetscErrorCode GetQuadratureDataP2C(PetscFE fe, PetscQuadrature quadrature, PetscInt *q_dim, PetscInt *num_comp, PetscInt *Q, 645 CeedScalar **q_points, CeedScalar **q_weights) { 646 PetscInt spatial_dim, dim, num_comp_quadrature, num_q_points; 647 const PetscReal *q_points_petsc, *q_weights_petsc; 648 649 PetscFunctionBeginUser; 650 PetscCall(PetscFEGetSpatialDimension(fe, &spatial_dim)); 651 PetscCall(PetscQuadratureGetData(quadrature, &dim, &num_comp_quadrature, &num_q_points, &q_points_petsc, &q_weights_petsc)); 652 if (q_dim) *q_dim = dim; 653 if (num_comp) *num_comp = num_comp_quadrature; 654 if (Q) *Q = num_q_points; 655 if (q_weights) { 656 PetscSizeT q_weights_size = num_q_points; 657 658 PetscCall(PetscCalloc1(q_weights_size, q_weights)); 659 for (CeedInt i = 0; i < num_q_points; i++) (*q_weights)[i] = q_weights_petsc[i]; 660 } 661 if (q_points) { 662 PetscSizeT q_points_size = num_q_points * spatial_dim; 663 664 PetscCall(PetscCalloc1(q_points_size, q_points)); 665 for (CeedInt i = 0; i < num_q_points; i++) { 666 for (CeedInt d = 0; d < dim; d++) (*q_points)[i * spatial_dim + d] = q_points_petsc[i * dim + d]; 667 } 668 } 669 PetscFunctionReturn(PETSC_SUCCESS); 670 } 671 672 /** 673 @brief Create 1D tabulation from `PetscFE`. 674 675 Collective across MPI processes. 676 677 @note `q_weights` and `q_points` are allocated using `PetscCalloc1` and must be freed by the user. 678 679 @param[in] comm `MPI_Comm` for creating `FE` object from options 680 @param[in] fe PETSc basis `FE` object 681 @param[out] tabulation PETSc basis tabulation 682 @param[out] q_points_1d Quadrature points in libCEED orientation 683 @param[out] q_weights_1d Quadrature weights in libCEED orientation 684 685 @return An error code: 0 - success, otherwise - failure 686 **/ 687 static inline PetscErrorCode Create1DTabulation_Tensor(MPI_Comm comm, PetscFE fe, PetscTabulation *tabulation, PetscReal **q_points_1d, 688 CeedScalar **q_weights_1d) { 689 PetscFE fe_1d; 690 PetscInt dim, num_comp, Q = -1, q_dim = -1, num_derivatives = 1; 691 PetscQuadrature quadrature; 692 693 PetscFunctionBeginUser; 694 // Create 1D FE 695 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 696 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 697 { 698 const char *prefix; 699 PetscBool is_tensor; 700 PetscInt num_comp, order, q_order; 701 PetscDualSpace dual_space; 702 PetscQuadrature quadrature; 703 704 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix)); 705 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 706 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 707 PetscCall(PetscDualSpaceLagrangeGetTensor(dual_space, &is_tensor)); 708 PetscCall(PetscDualSpaceGetOrder(dual_space, &order)); 709 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 710 PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &q_order, NULL, NULL)); 711 PetscCall(PetscFECreateLagrangeFromOptions(comm, 1, num_comp, !is_tensor, order, 712 (PetscInt)PetscCeilReal(PetscPowReal(1.0 * q_order, 1.0 / dim)) - 1, prefix, &fe_1d)); 713 } 714 715 // Get quadrature data 716 PetscCall(PetscFEGetQuadrature(fe_1d, &quadrature)); 717 PetscCall(GetQuadratureDataP2C(fe_1d, quadrature, &q_dim, NULL, &Q, q_points_1d, q_weights_1d)); 718 719 // Compute 1D tabulation 720 PetscCall(PetscFECreateTabulation(fe_1d, 1, Q, *q_points_1d, num_derivatives, tabulation)); 721 722 // Cleanup 723 PetscCall(PetscFEDestroy(&fe_1d)); 724 725 PetscFunctionReturn(PETSC_SUCCESS); 726 } 727 728 /** 729 @brief Destroy `CeedBasis` in a `PetscContainer`. 730 731 Not collective across MPI processes. 732 733 @param[in,out] ctx `CeedBasis` 734 735 @return An error code: 0 - success, otherwise - failure 736 **/ 737 static PetscErrorCode DMPlexCeedBasisDestroy(void **ctx) { 738 CeedBasis basis = *(CeedBasis *)ctx; 739 740 PetscFunctionBegin; 741 PetscCallCeed(CeedBasisReturnCeed(basis), CeedBasisDestroy(&basis)); 742 PetscFunctionReturn(PETSC_SUCCESS); 743 } 744 745 /** 746 @brief Create `CeedBasis` from `DMPlex`. 747 748 The `CeedBasis` is stored in the `DMPlex` object for reuse; 749 if the routine is called again with the same arguments, the same `CeedBasis` is returned. 750 751 Collective across MPI processes. 752 753 @param[in] ceed `Ceed` context 754 @param[in] dm `DMPlex` holding mesh 755 @param[in] domain_label `DMLabel` for `DMPlex` domain 756 @param[in] label_value Stratum value 757 @param[in] height Height of `DMPlex` topology 758 @param[in] dm_field Index of `DMPlex` field 759 @param[out] basis `CeedBasis` for `DMPlex` 760 761 @return An error code: 0 - success, otherwise - failure 762 **/ 763 PetscErrorCode DMPlexCeedBasisCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 764 CeedBasis *basis) { 765 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 766 char container_name[1024]; 767 CeedBasis container_basis; 768 769 PetscFunctionBeginUser; 770 { 771 const char *label_name = NULL; 772 773 if (domain_label) PetscCall(PetscObjectGetName((PetscObject)domain_label, &label_name)); 774 PetscCall(PetscSNPrintf(container_name, sizeof(container_name), 775 "DM CeedBasis - label: %s label value: %" PetscInt_FMT " height: %" PetscInt_FMT " DM field: %" PetscInt_FMT, 776 label_name ? label_name : "NONE", label_value, height, dm_field)); 777 } 778 PetscCall(PetscObjectContainerQuery((PetscObject)dm, container_name, (void **)&container_basis)); 779 780 if (container_basis) { 781 // Copy existing basis 782 *basis = NULL; 783 PetscCallCeed(ceed, CeedBasisReferenceCopy(container_basis, basis)); 784 } else { 785 PetscDS ds; 786 PetscFE fe; 787 PetscDualSpace dual_space; 788 PetscQuadrature quadrature; 789 PetscBool is_tensor = PETSC_TRUE; 790 PetscInt ds_field = -1; 791 792 // Get element information 793 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 794 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 795 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 796 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 797 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 798 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 799 PetscCall(PetscDualSpaceLagrangeGetTensor(dual_space, &is_tensor)); 800 801 // Build libCEED basis 802 PetscBool use_nontensor = PETSC_FALSE; 803 804 PetscOptionsBegin(comm, NULL, "Tensor Basis as Nontensor Check", NULL); 805 PetscCall(PetscOptionsBool("-ceed_basis_as_nontensor", "", NULL, use_nontensor, &use_nontensor, NULL)); 806 PetscOptionsEnd(); 807 808 if (!is_tensor || use_nontensor) { 809 PetscTabulation basis_tabulation; 810 PetscInt num_derivatives = 1, face = 0; 811 CeedScalar *q_points, *q_weights, *interp, *grad; 812 CeedElemTopology elem_topo; 813 814 { 815 DMPolytopeType cell_type; 816 817 PetscCall(GetGlobalDMPlexPolytopeType(dm, domain_label, label_value, height, &cell_type)); 818 elem_topo = PolytopeTypePetscToCeed(cell_type); 819 PetscCheck(elem_topo, comm, PETSC_ERR_SUP, "DMPlex topology not supported: %s", DMPolytopeTypes[cell_type]); 820 } 821 822 // Compute basis tabulation 823 PetscCall(GetQuadratureDataP2C(fe, quadrature, NULL, NULL, NULL, &q_points, &q_weights)); 824 PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation)); 825 PetscCall(ComputeFieldTabulationP2C(dm, dm_field, face, basis_tabulation, &interp, &grad)); 826 { 827 PetscInt num_comp = basis_tabulation->Nc, P = basis_tabulation->Nb / num_comp, Q = basis_tabulation->Np; 828 829 PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis)); 830 } 831 832 PetscCall(PetscFree(q_points)); 833 PetscCall(PetscFree(q_weights)); 834 PetscCall(PetscFree(interp)); 835 PetscCall(PetscFree(grad)); 836 } else { 837 PetscInt P_1d, Q_1d, num_comp, dim; 838 PetscTabulation basis_tabulation; 839 CeedScalar *q_points_1d, *q_weights_1d, *interp_1d, *grad_1d; 840 841 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 842 PetscCall(Create1DTabulation_Tensor(comm, fe, &basis_tabulation, &q_points_1d, &q_weights_1d)); 843 num_comp = basis_tabulation->Nc; 844 P_1d = basis_tabulation->Nb / num_comp; 845 Q_1d = basis_tabulation->Np; 846 PetscCall(ComputeFieldTabulationP2C(dm, dm_field, 0, basis_tabulation, &interp_1d, &grad_1d)); 847 PetscCallCeed(ceed, CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_points_1d, q_weights_1d, basis)); 848 849 // Cleanup 850 PetscCall(PetscFree(q_points_1d)); 851 PetscCall(PetscFree(q_weights_1d)); 852 PetscCall(PetscFree(interp_1d)); 853 PetscCall(PetscFree(grad_1d)); 854 PetscCall(PetscTabulationDestroy(&basis_tabulation)); 855 } 856 // Set in container 857 PetscCallCeed(ceed, CeedBasisReferenceCopy(*basis, &container_basis)); 858 PetscCall(PetscObjectContainerCompose((PetscObject)dm, container_name, container_basis, DMPlexCeedBasisDestroy)); 859 } 860 PetscFunctionReturn(PETSC_SUCCESS); 861 } 862 863 /** 864 @brief Create `CeedBasis` for coordinate space of `DMPlex`. 865 866 The `CeedBasis` is stored in the coordinate `DMPlex` object for reuse; 867 if the routine is called again with the same arguments, the same `CeedBasis` is returned. 868 869 Collective across MPI processes. 870 871 @param[in] ceed `Ceed` context 872 @param[in] dm `DMPlex` holding mesh 873 @param[in] domain_label `DMLabel` for `DMPlex` domain 874 @param[in] label_value Stratum value 875 @param[in] height Height of `DMPlex` topology 876 @param[out] basis `CeedBasis` for coordinate space of `DMPlex` 877 **/ 878 PetscErrorCode DMPlexCeedBasisCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, CeedBasis *basis) { 879 DM dm_coord = NULL; 880 881 PetscFunctionBeginUser; 882 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 883 PetscCall(DMPlexCeedBasisCreate(ceed, dm_coord, domain_label, label_value, height, 0, basis)); 884 PetscFunctionReturn(PETSC_SUCCESS); 885 } 886 887 typedef struct CeedVectorAndState_ *CeedVectorAndState; 888 struct CeedVectorAndState_ { 889 CeedVector vector; // !< `CeedVector` containing copy of data from a `Vec` 890 PetscObjectState state; // !< Object state of the `Vec` when the data was copied 891 }; 892 893 /** 894 @brief Destroy `CeedVectorAndState` in a `PetscContainer`. 895 896 Not collective across MPI processes. 897 898 @param[in,out] ctx `CeedVectorAndState` 899 900 @return An error code: 0 - success, otherwise - failure 901 **/ 902 static PetscErrorCode CeedVectorAndStateDestroy(void **ctx) { 903 CeedVectorAndState vec_and_state = *(CeedVectorAndState *)ctx; 904 905 PetscFunctionBegin; 906 PetscCallCeed(CeedVectorReturnCeed(vec_and_state->vector), CeedVectorDestroy(&vec_and_state->vector)); 907 PetscCall(PetscFree(vec_and_state)); 908 *ctx = NULL; 909 PetscFunctionReturn(PETSC_SUCCESS); 910 } 911 912 /** 913 @brief Get Ceed objects for coordinate field of DM 914 915 @param[in] ceed `Ceed` context 916 @param[in] dm `DMPlex` holding mesh 917 @param[in] domain_label `DMLabel` for `DMPlex` domain 918 @param[in] label_value Stratum value 919 @param[in] height Height of `DMPlex` topology 920 @param[out] elem_restr `CeedElemRestriction` for coodinates of `DMPlex`, or `NULL` 921 @param[out] basis `CeedBasis` for coordinate space of `DMPlex`, or `NULL` 922 @param[out] vector `CeedVector` for coordinates of `DMPlex`, or `NULL` 923 **/ 924 PetscErrorCode DMPlexCeedCoordinateCreateField(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 925 CeedElemRestriction *elem_restr, CeedBasis *basis, CeedVector *vector) { 926 CeedElemRestriction elem_restr_ = NULL; 927 928 PetscFunctionBeginUser; 929 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_)); 930 if (elem_restr) { 931 *elem_restr = NULL; 932 PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(elem_restr_, elem_restr)); 933 } 934 if (basis) PetscCall(DMPlexCeedBasisCoordinateCreate(ceed, dm, domain_label, label_value, height, basis)); 935 if (vector) { 936 DM cdm; 937 char container_name[] = "DM Coordinate CeedVectorAndState"; 938 CeedVectorAndState vec_and_state; 939 Vec X_loc; 940 PetscObjectState X_loc_state; 941 942 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 943 if (cdm) { 944 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 945 } else { 946 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 947 } 948 PetscCall(VecGetState(X_loc, &X_loc_state)); 949 PetscCall(PetscObjectContainerQuery((PetscObject)X_loc, container_name, (void **)&vec_and_state)); 950 951 if (vec_and_state && vec_and_state->state == X_loc_state) { 952 // Copy existing vector 953 *vector = NULL; 954 PetscCallCeed(ceed, CeedVectorReferenceCopy(vec_and_state->vector, vector)); 955 } else { 956 PetscCall(CeedElemRestrictionCreateVector(elem_restr_, vector, NULL)); 957 PetscCall(VecCopyPetscToCeed(X_loc, *vector)); 958 959 // Set in container 960 PetscCall(PetscNew(&vec_and_state)); 961 PetscCallCeed(ceed, CeedVectorReferenceCopy(*vector, &vec_and_state->vector)); 962 vec_and_state->state = X_loc_state; 963 PetscCall(PetscObjectContainerCompose((PetscObject)X_loc, container_name, vec_and_state, CeedVectorAndStateDestroy)); 964 } 965 } 966 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_)); 967 PetscFunctionReturn(PETSC_SUCCESS); 968 } 969 970 /** 971 @brief Create `CeedBasis` for cell to face quadrature space evaluation from `DMPlex`. 972 973 Collective across MPI processes. 974 975 @param[in] ceed `Ceed` context 976 @param[in] dm `DMPlex` holding mesh 977 @param[in] domain_label `DMLabel` for `DMPlex` domain 978 @param[in] label_value Stratum value 979 @param[in] face Index of face 980 @param[in] dm_field Index of `DMPlex` field 981 @param[out] basis `CeedBasis` for cell to face evaluation 982 **/ 983 PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field, 984 CeedBasis *basis) { 985 PetscInt ds_field = -1, height = 0; 986 PetscDS ds; 987 PetscFE fe; 988 PetscQuadrature face_quadrature; 989 990 PetscFunctionBeginUser; 991 // Get element information 992 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 993 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 994 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 995 PetscCall(PetscFEGetFaceQuadrature(fe, &face_quadrature)); 996 997 { // Build libCEED basis 998 PetscInt num_derivatives = 1, num_comp, P, Q = -1; 999 PetscTabulation basis_tabulation; 1000 CeedScalar *q_points, *q_weights, *interp, *grad; 1001 CeedElemTopology elem_topo; 1002 1003 { // Verify compatible element topology 1004 DMPolytopeType cell_type; 1005 1006 PetscCall(GetGlobalDMPlexPolytopeType(dm, domain_label, label_value, height, &cell_type)); 1007 elem_topo = PolytopeTypePetscToCeed(cell_type); 1008 PetscCheck(elem_topo, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported: %s", DMPolytopeTypes[cell_type]); 1009 } 1010 1011 // Compute basis tabulation 1012 PetscCall(GetQuadratureDataP2C(fe, face_quadrature, NULL, NULL, &Q, &q_points, &q_weights)); 1013 PetscCall(PetscFEGetFaceTabulation(fe, num_derivatives, &basis_tabulation)); 1014 num_comp = basis_tabulation->Nc; 1015 P = basis_tabulation->Nb / num_comp; 1016 PetscCall(ComputeFieldTabulationP2C(dm, dm_field, face, basis_tabulation, &interp, &grad)); 1017 PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis)); 1018 1019 PetscCall(PetscFree(q_points)); 1020 PetscCall(PetscFree(q_weights)); 1021 PetscCall(PetscFree(interp)); 1022 PetscCall(PetscFree(grad)); 1023 } 1024 PetscFunctionReturn(PETSC_SUCCESS); 1025 } 1026 1027 /** 1028 @brief Create `CeedBasis` for cell to face quadrature space evaluation on coordinate space of `DMPlex`. 1029 1030 Collective across MPI processes. 1031 1032 @param[in] ceed `Ceed` context 1033 @param[in] dm `DMPlex` holding mesh 1034 @param[in] domain_label `DMLabel` for `DMPlex` domain 1035 @param[in] label_value Stratum value 1036 @param[in] face Index of face 1037 @param[out] basis `CeedBasis` for cell to face evaluation on coordinate space of `DMPlex` 1038 **/ 1039 PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, 1040 CeedBasis *basis) { 1041 DM dm_coord; 1042 1043 PetscFunctionBeginUser; 1044 PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 1045 if (!dm_coord) { 1046 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1047 } 1048 PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, dm_coord, domain_label, label_value, face, 0, basis)); 1049 PetscFunctionReturn(PETSC_SUCCESS); 1050 } 1051 1052 /** 1053 @brief Setup `DM` with FE space of appropriate degree 1054 1055 Must be followed by `DMSetupByOrderEnd_FEM` 1056 1057 @param[in] setup_faces Flag to setup face geometry 1058 @param[in] setup_coords Flag to setup coordinate spaces 1059 @param[in] degree Polynomial orders of field 1060 @param[in] coord_order Polynomial order of coordinate basis, or `PETSC_DECIDE` for default 1061 @param[in] q_extra Additional quadrature order 1062 @param[in] num_fields Number of fields in solution vector 1063 @param[in] field_sizes Array of number of components for each field 1064 @param[out] dm `DM` to setup 1065 1066 @return An error code: 0 - success, otherwise - failure 1067 **/ 1068 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 1069 PetscInt num_fields, const PetscInt *field_sizes, DM dm) { 1070 PetscInt dim, q_order = degree + q_extra; 1071 PetscBool is_simplex = PETSC_TRUE; 1072 PetscFE fe; 1073 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 1074 1075 PetscFunctionBeginUser; 1076 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 1077 1078 // Setup DM 1079 PetscCall(DMGetDimension(dm, &dim)); 1080 for (PetscInt i = 0; i < num_fields; i++) { 1081 PetscFE fe_face; 1082 PetscInt q_order = degree + q_extra; 1083 1084 PetscCall(PetscFECreateLagrange(comm, dim, field_sizes[i], is_simplex, degree, q_order, &fe)); 1085 if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe_face)); 1086 PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 1087 PetscCall(PetscFEDestroy(&fe)); 1088 } 1089 PetscCall(DMCreateDS(dm)); 1090 1091 // Project coordinates to enrich quadrature space 1092 if (setup_coords) { 1093 DM dm_coord; 1094 PetscDS ds_coord; 1095 PetscFE fe_coord_current, fe_coord_new, fe_coord_face_new; 1096 PetscDualSpace fe_coord_dual_space; 1097 PetscInt fe_coord_order, num_comp_coord; 1098 1099 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1100 PetscCall(DMGetCoordinateDim(dm, &num_comp_coord)); 1101 PetscCall(DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord, NULL)); 1102 PetscCall(PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current)); 1103 PetscCall(PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space)); 1104 PetscCall(PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order)); 1105 1106 // Create FE for coordinates 1107 if (coord_order != PETSC_DECIDE) fe_coord_order = coord_order; 1108 PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, q_order, &fe_coord_new)); 1109 if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe_coord_new, 1, &fe_coord_face_new)); 1110 PetscCall(DMSetCoordinateDisc(dm, fe_coord_new, PETSC_FALSE, PETSC_TRUE)); 1111 PetscCall(PetscFEDestroy(&fe_coord_new)); 1112 } 1113 PetscFunctionReturn(PETSC_SUCCESS); 1114 } 1115 1116 /** 1117 @brief Finish setting up `DM` 1118 1119 Must be called after `DMSetupByOrderBegin_FEM` and all strong BCs have be declared via `DMAddBoundaries`. 1120 1121 @param[in] setup_coords Flag to setup coordinate spaces 1122 @param[out] dm `DM` to setup 1123 1124 @return An error code: 0 - success, otherwise - failure 1125 **/ 1126 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm) { 1127 PetscBool is_simplex; 1128 1129 PetscFunctionBeginUser; 1130 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 1131 // Set tensor permutation if needed 1132 if (!is_simplex) { 1133 PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 1134 if (setup_coords) { 1135 DM dm_coord; 1136 1137 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1138 PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 1139 } 1140 } 1141 PetscCall(DMLocalizeCoordinates(dm)); // Must localize after tensor closure setting 1142 PetscFunctionReturn(PETSC_SUCCESS); 1143 } 1144 1145 /** 1146 @brief Setup `DM` with FE space of appropriate degree with no boundary conditions 1147 1148 Calls `DMSetupByOrderBegin_FEM` and `DMSetupByOrderEnd_FEM` successively 1149 1150 @param[in] setup_faces Flag to setup face geometry 1151 @param[in] setup_coords Flag to setup coordinate spaces 1152 @param[in] degree Polynomial orders of field 1153 @param[in] coord_order Polynomial order of coordinate basis, or `PETSC_DECIDE` for default 1154 @param[in] q_extra Additional quadrature order 1155 @param[in] num_fields Number of fields in solution vector 1156 @param[in] field_sizes Array of number of components for each field 1157 @param[out] dm `DM` to setup 1158 1159 @return An error code: 0 - success, otherwise - failure 1160 **/ 1161 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 1162 PetscInt num_fields, const PetscInt *field_sizes, DM dm) { 1163 PetscFunctionBeginUser; 1164 PetscCall(DMSetupByOrderBegin_FEM(setup_faces, setup_coords, degree, coord_order, q_extra, num_fields, field_sizes, dm)); 1165 PetscCall(DMSetupByOrderEnd_FEM(setup_coords, dm)); 1166 PetscFunctionReturn(PETSC_SUCCESS); 1167 } 1168 1169 /** 1170 @brief Get the points in a label stratum which have requested height 1171 1172 Example, from the "Face Sets" label (which contains face, edge, and vertex points), return points that correspond to face points (height = 1) and have a label value of 3: 1173 ``` 1174 PetscCall(DMGetStratumISAtHeight(dm, "Face Sets", 3, 1, &is_face_points)); 1175 ``` 1176 1177 @param[in] dm DM which has the label 1178 @param[in] name Name fo the label 1179 @param[in] value Label value of points to return 1180 @param[in] height Height of points to return 1181 @param[out] points `IS` of points, or `NULL` is there are no points in that stratum at height 1182 **/ 1183 static PetscErrorCode DMGetStratumISAtHeight(DM dm, const char name[], PetscInt value, PetscInt height, IS *points) { 1184 PetscInt depth; 1185 DMLabel depth_label; 1186 IS label_is, depth_is; 1187 1188 PetscFunctionBeginUser; 1189 PetscCall(DMPlexGetDepth(dm, &depth)); 1190 PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 1191 PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 1192 PetscCall(DMGetStratumIS(dm, name, value, &label_is)); 1193 if (label_is == NULL || depth_is == NULL) *points = NULL; 1194 else PetscCall(ISIntersect(depth_is, label_is, points)); 1195 PetscCall(ISDestroy(&depth_is)); 1196 PetscCall(ISDestroy(&label_is)); 1197 PetscFunctionReturn(PETSC_SUCCESS); 1198 } 1199 1200 /** 1201 @brief Create a label with separate values for the `FE` orientations for all cells with this material for the `DMPlex` face. 1202 1203 Collective across MPI processes. 1204 1205 @param[in,out] dm `DMPlex` holding mesh 1206 @param[in] dm_face Face number on `DMPlex` 1207 @param[out] face_label_name Label name for newly created label, or `NULL` if no elements match the `material` and `dm_face` 1208 **/ 1209 PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, char **face_label_name) { 1210 PetscInt num_face_points = 0, face_height = 1; 1211 const PetscInt *face_points; 1212 IS is_face_points; 1213 DMLabel face_label = NULL; 1214 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 1215 1216 PetscFunctionBeginUser; 1217 { // Create label name and label 1218 PetscSizeT label_name_len = 64; 1219 PetscBool has_label_already; 1220 1221 // Face label 1222 PetscCall(PetscCalloc1(label_name_len, face_label_name)); 1223 PetscCall(PetscSNPrintf(*face_label_name, label_name_len, "Orientation for DM face %" PetscInt_FMT, dm_face)); 1224 PetscCall(DMHasLabel(dm, *face_label_name, &has_label_already)); 1225 if (has_label_already) PetscFunctionReturn(PETSC_SUCCESS); 1226 1227 PetscCall(DMCreateLabel(dm, *face_label_name)); 1228 PetscCall(DMGetLabel(dm, *face_label_name, &face_label)); 1229 } 1230 1231 PetscCall(DMGetStratumISAtHeight(dm, "Face Sets", dm_face, face_height, &is_face_points)); 1232 if (is_face_points) { 1233 PetscCall(ISGetLocalSize(is_face_points, &num_face_points)); 1234 PetscCall(ISGetIndices(is_face_points, &face_points)); 1235 } 1236 1237 for (PetscInt face_index = 0; face_index < num_face_points; face_index++) { 1238 const PetscInt *face_support, *cell_cone; 1239 PetscInt face_support_size = 0, cell_cone_size = 0, fe_face_on_cell = -1; 1240 PetscInt face_point = face_points[face_index]; 1241 1242 // Get cell supporting face 1243 PetscCall(DMPlexGetSupport(dm, face_point, &face_support)); 1244 PetscCall(DMPlexGetSupportSize(dm, face_point, &face_support_size)); 1245 PetscCheck(face_support_size == 1, comm, PETSC_ERR_SUP, "Expected one cell in support of exterior face, but got %" PetscInt_FMT " cells", 1246 face_support_size); 1247 PetscInt cell_point = face_support[0]; 1248 1249 // Find face number 1250 PetscCall(DMPlexGetCone(dm, cell_point, &cell_cone)); 1251 PetscCall(DMPlexGetConeSize(dm, cell_point, &cell_cone_size)); 1252 for (PetscInt i = 0; i < cell_cone_size; i++) { 1253 if (cell_cone[i] == face_point) fe_face_on_cell = i; 1254 } 1255 PetscCheck(fe_face_on_cell != -1, comm, PETSC_ERR_SUP, "Could not find face %" PetscInt_FMT " in cone of its support", face_point); 1256 1257 // Add to label 1258 PetscCall(DMLabelSetValue(face_label, face_point, fe_face_on_cell)); 1259 } 1260 PetscCall(DMPlexLabelAddFaceCells(dm, face_label)); 1261 PetscCall(ISDestroy(&is_face_points)); 1262 PetscFunctionReturn(PETSC_SUCCESS); 1263 } 1264 1265 /** 1266 @brief Creates an array of all the values set for a `DMLabel` 1267 1268 Because `DMLabelGetValueIS` returns label values locally, not globally. 1269 1270 Collective across MPI processes. 1271 1272 @param[in] dm `DM` with the label 1273 @param[in] label `DMLabel` to get values of 1274 @param[out] num_values Total number of values 1275 @param[out] value_array Array of label values, must be freed by user 1276 **/ 1277 PetscErrorCode DMLabelCreateGlobalValueArray(DM dm, DMLabel label, PetscInt *num_values, PetscInt **value_array) { 1278 PetscInt num_values_local, minmax_values[2], minmax_values_loc[2]; 1279 PetscInt *values_local = NULL; 1280 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 1281 PetscSegBuffer seg; 1282 PetscCount seg_size; 1283 1284 PetscFunctionBeginUser; 1285 { // Extract array of local DMLabel values so it can be sorted 1286 IS is_values; 1287 const PetscInt *is_values_local = NULL; 1288 1289 PetscCall(DMLabelGetValueIS(label, &is_values)); 1290 PetscCall(ISGetLocalSize(is_values, &num_values_local)); 1291 1292 PetscCall(ISGetIndices(is_values, &is_values_local)); 1293 PetscCall(PetscMalloc1(num_values_local, &values_local)); 1294 PetscCall(PetscArraycpy(values_local, is_values_local, num_values_local)); 1295 PetscCall(ISRestoreIndices(is_values, &is_values_local)); 1296 PetscCall(ISDestroy(&is_values)); 1297 } 1298 1299 if (num_values_local) { 1300 PetscCall(PetscSortInt(num_values_local, values_local)); 1301 minmax_values_loc[0] = values_local[0]; 1302 minmax_values_loc[1] = values_local[num_values_local - 1]; 1303 } else { 1304 // Rank has no set label values, therefore set local min/max to have no effect on global min/max 1305 minmax_values_loc[0] = PETSC_INT_MAX; 1306 minmax_values_loc[1] = PETSC_INT_MIN; 1307 } 1308 PetscCall(PetscGlobalMinMaxInt(comm, minmax_values_loc, minmax_values)); 1309 1310 PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 16, &seg)); 1311 for (PetscInt value = minmax_values[0]; value <= minmax_values[1]; value++) { 1312 PetscInt location_local = -1, location_global = -1; 1313 1314 PetscCall(PetscFindInt(value, num_values_local, values_local, &location_local)); 1315 PetscCallMPI(MPIU_Allreduce(&location_local, &location_global, 1, MPIU_INT, MPI_MAX, comm)); 1316 if (location_global < 0) continue; 1317 else { 1318 PetscInt *seg_slot; 1319 1320 PetscCall(PetscSegBufferGetInts(seg, 1, &seg_slot)); 1321 *seg_slot = value; 1322 } 1323 } 1324 PetscCall(PetscSegBufferGetSize(seg, &seg_size)); 1325 *num_values = (PetscInt)seg_size; 1326 PetscCall(PetscSegBufferExtractAlloc(seg, value_array)); 1327 PetscCall(PetscSegBufferDestroy(&seg)); 1328 if (values_local) PetscCall(PetscFree(values_local)); 1329 PetscFunctionReturn(PETSC_SUCCESS); 1330 } 1331 1332 /** 1333 @brief Get the number of components for `dm`'s coordinate field 1334 1335 @param[in] dm DM 1336 @param[out] num_comp Number of components for the coordinate field 1337 **/ 1338 PetscErrorCode DMGetCoordinateNumComps(DM dm, PetscInt *num_comp) { 1339 DM dm_coord; 1340 PetscSection section_coord; 1341 PetscInt field = 0; // Default field has the coordinates 1342 1343 PetscFunctionBeginUser; 1344 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1345 PetscCall(DMGetLocalSection(dm_coord, §ion_coord)); 1346 PetscCall(PetscSectionGetFieldComponents(section_coord, field, num_comp)); 1347 PetscFunctionReturn(PETSC_SUCCESS); 1348 } 1349 1350 /** 1351 @brief Get the number of components for `dm`'s field 1352 1353 @param[in] dm DM 1354 @param[in] field The field ID to get the number of components for 1355 @param[out] num_comp Number of components for the coordinate field 1356 **/ 1357 PetscErrorCode DMGetFieldNumComps(DM dm, PetscInt field, PetscInt *num_comp) { 1358 PetscSection section; 1359 1360 PetscFunctionBeginUser; 1361 PetscCall(DMGetLocalSection(dm, §ion)); 1362 PetscCall(PetscSectionGetFieldComponents(section, field, num_comp)); 1363 PetscFunctionReturn(PETSC_SUCCESS); 1364 } 1365