1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utilities for setting up DM and PetscFE 10 11 #include <ceed.h> 12 #include <petscdmplex.h> 13 #include <petscds.h> 14 15 #include "../navierstokes.h" 16 17 /** 18 @brief Convert `DM` field to `DS` field. 19 20 @param[in] dm `DM` holding mesh 21 @param[in] domain_label Label for `DM` domain 22 @param[in] dm_field Index of `DM` field 23 @param[out] ds_field Index of `DS` field 24 25 @return An error code: 0 - success, otherwise - failure 26 **/ 27 PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) { 28 PetscDS ds; 29 IS field_is; 30 const PetscInt *fields; 31 PetscInt num_fields; 32 33 PetscFunctionBeginUser; 34 PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL)); 35 PetscCall(ISGetIndices(field_is, &fields)); 36 PetscCall(ISGetSize(field_is, &num_fields)); 37 for (PetscInt i = 0; i < num_fields; i++) { 38 if (dm_field == fields[i]) { 39 *ds_field = i; 40 break; 41 } 42 } 43 PetscCall(ISRestoreIndices(field_is, &fields)); 44 45 if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 46 PetscFunctionReturn(PETSC_SUCCESS); 47 } 48 49 /** 50 @brief Create `CeedElemRestriction` from `DMPlex`. 51 52 Not collective across MPI processes. 53 54 @param[in] ceed `Ceed` context 55 @param[in] dm `DMPlex` holding mesh 56 @param[in] domain_label `DMLabel` for `DMPlex` domain 57 @param[in] label_value Stratum value 58 @param[in] height Height of `DMPlex` topology 59 @param[in] dm_field Index of `DMPlex` field 60 @param[out] restriction `CeedElemRestriction` for `DMPlex` 61 62 @return An error code: 0 - success, otherwise - failure 63 **/ 64 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 65 CeedElemRestriction *restriction) { 66 PetscInt num_elem, elem_size, num_dof, num_comp, *restriction_offsets_petsc; 67 CeedInt *restriction_offsets_ceed = NULL; 68 69 PetscFunctionBeginUser; 70 PetscCall( 71 DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &restriction_offsets_petsc)); 72 PetscCall(IntArrayP2C(num_elem * elem_size, &restriction_offsets_petsc, &restriction_offsets_ceed)); 73 PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 74 restriction_offsets_ceed, restriction)); 75 PetscCall(PetscFree(restriction_offsets_ceed)); 76 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79 80 /** 81 @brief Create `CeedElemRestriction` from `DMPlex` domain for mesh coordinates. 82 83 Not collective across MPI processes. 84 85 @param[in] ceed `Ceed` context 86 @param[in] dm `DMPlex` holding mesh 87 @param[in] domain_label Label for `DMPlex` domain 88 @param[in] label_value Stratum value 89 @param[in] height Height of `DMPlex` topology 90 @param[out] restriction `CeedElemRestriction` for mesh 91 92 @return An error code: 0 - success, otherwise - failure 93 **/ 94 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 95 CeedElemRestriction *restriction) { 96 DM dm_coord; 97 PetscBool is_simplex = PETSC_FALSE; 98 99 PetscFunctionBeginUser; 100 PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 101 if (!dm_coord) { 102 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 103 } 104 PetscCall(DMPlexIsSimplex(dm_coord, &is_simplex)); 105 if (!is_simplex) PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 106 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm_coord, domain_label, label_value, height, 0, restriction)); 107 108 PetscFunctionReturn(PETSC_SUCCESS); 109 } 110 111 /** 112 @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data. 113 114 Not collective across MPI processes. 115 116 @param[in] ceed `Ceed` context 117 @param[in] dm `DMPlex` holding mesh 118 @param[in] domain_label Label for `DMPlex` domain 119 @param[in] label_value Stratum value 120 @param[in] height Height of `DMPlex` topology 121 @param[in] q_data_size Number of components for `QFunction` data 122 @param[in] is_collocated Boolean flag indicating if the data is collocated on the nodes (`PETSC_TRUE`) on on quadrature points (`PETSC_FALSE`) 123 @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 124 125 @return An error code: 0 - success, otherwise - failure 126 **/ 127 static PetscErrorCode DMPlexCeedElemRestrictionStridedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 128 PetscInt q_data_size, PetscBool is_collocated, CeedElemRestriction *restriction) { 129 PetscInt num_elem, num_qpts, dm_field = 0; 130 131 PetscFunctionBeginUser; 132 { // Get number of elements 133 PetscInt depth; 134 DMLabel depth_label; 135 IS point_is, depth_is; 136 137 PetscCall(DMPlexGetDepth(dm, &depth)); 138 PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 139 PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 140 if (domain_label) { 141 IS domain_is; 142 143 PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is)); 144 if (domain_is) { 145 PetscCall(ISIntersect(depth_is, domain_is, &point_is)); 146 PetscCall(ISDestroy(&domain_is)); 147 } else { 148 point_is = NULL; 149 } 150 PetscCall(ISDestroy(&depth_is)); 151 } else { 152 point_is = depth_is; 153 } 154 if (point_is) { 155 PetscCall(ISGetLocalSize(point_is, &num_elem)); 156 } else { 157 num_elem = 0; 158 } 159 PetscCall(ISDestroy(&point_is)); 160 } 161 162 { // Get number of quadrature points 163 PetscDS ds; 164 PetscFE fe; 165 PetscInt ds_field = -1; 166 167 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 168 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 169 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 170 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 171 if (is_collocated) { 172 PetscDualSpace dual_space; 173 PetscInt num_dual_basis_vectors, dim, num_comp; 174 175 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 176 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 177 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 178 PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 179 num_qpts = num_dual_basis_vectors / num_comp; 180 } else { 181 PetscQuadrature quadrature; 182 183 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 184 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 185 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 186 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 187 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 188 PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &num_qpts, NULL, NULL)); 189 } 190 } 191 192 // Create the restriction 193 PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, num_elem * num_qpts * q_data_size, CEED_STRIDES_BACKEND, 194 restriction)); 195 196 PetscFunctionReturn(PETSC_SUCCESS); 197 } 198 199 /** 200 @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data. 201 202 Not collective across MPI processes. 203 204 @param[in] ceed `Ceed` context 205 @param[in] dm `DMPlex` holding mesh 206 @param[in] domain_label Label for `DMPlex` domain 207 @param[in] label_value Stratum value 208 @param[in] height Height of `DMPlex` topology 209 @param[in] q_data_size Number of components for `QFunction` data 210 @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 211 212 @return An error code: 0 - success, otherwise - failure 213 **/ 214 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 215 PetscInt q_data_size, CeedElemRestriction *restriction) { 216 PetscFunctionBeginUser; 217 PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_FALSE, restriction)); 218 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221 222 /** 223 @brief Create `CeedElemRestriction` from `DMPlex` domain for nodally collocated auxilury `QFunction` data. 224 225 Not collective across MPI processes. 226 227 @param[in] ceed `Ceed` context 228 @param[in] dm `DMPlex` holding mesh 229 @param[in] domain_label Label for `DMPlex` domain 230 @param[in] label_value Stratum value 231 @param[in] height Height of `DMPlex` topology 232 @param[in] q_data_size Number of components for `QFunction` data 233 @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 234 235 @return An error code: 0 - success, otherwise - failure 236 **/ 237 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 238 PetscInt q_data_size, CeedElemRestriction *restriction) { 239 PetscFunctionBeginUser; 240 PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_TRUE, restriction)); 241 242 PetscFunctionReturn(PETSC_SUCCESS); 243 } 244 245 // ----------------------------------------------------------------------------- 246 // Utility function - convert from DMPolytopeType to CeedElemTopology 247 // ----------------------------------------------------------------------------- 248 static CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) { 249 switch (cell_type) { 250 case DM_POLYTOPE_TRIANGLE: 251 return CEED_TOPOLOGY_TRIANGLE; 252 case DM_POLYTOPE_QUADRILATERAL: 253 return CEED_TOPOLOGY_QUAD; 254 case DM_POLYTOPE_TETRAHEDRON: 255 return CEED_TOPOLOGY_TET; 256 case DM_POLYTOPE_HEXAHEDRON: 257 return CEED_TOPOLOGY_HEX; 258 default: 259 return 0; 260 } 261 } 262 263 // ----------------------------------------------------------------------------- 264 // Create libCEED Basis from PetscTabulation 265 // ----------------------------------------------------------------------------- 266 PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe, 267 PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) { 268 PetscInt first_point; 269 PetscInt ids[1] = {label_value}; 270 DMLabel depth_label; 271 DMPolytopeType cell_type; 272 CeedElemTopology elem_topo; 273 PetscScalar *q_points, *interp, *grad; 274 const PetscScalar *q_weights; 275 PetscDualSpace dual_space; 276 PetscInt num_dual_basis_vectors; 277 PetscInt dim, num_comp, P, Q; 278 279 PetscFunctionBeginUser; 280 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 281 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 282 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 283 PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 284 P = num_dual_basis_vectors / num_comp; 285 286 // Use depth label if no domain label present 287 if (!domain_label) { 288 PetscInt depth; 289 290 PetscCall(DMPlexGetDepth(dm, &depth)); 291 PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 292 ids[0] = depth - height; 293 } 294 295 // Get cell interp, grad, and quadrature data 296 PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL)); 297 PetscCall(DMPlexGetCellType(dm, first_point, &cell_type)); 298 elem_topo = ElemTopologyP2C(cell_type); 299 if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported"); 300 { 301 size_t q_points_size; 302 const PetscScalar *q_points_petsc; 303 PetscInt q_dim; 304 305 PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights)); 306 q_points_size = Q * dim * sizeof(CeedScalar); 307 PetscCall(PetscCalloc(q_points_size, &q_points)); 308 for (PetscInt q = 0; q < Q; q++) { 309 for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d]; 310 } 311 } 312 313 { // Convert to libCEED orientation 314 PetscBool is_simplex = PETSC_FALSE; 315 IS permutation = NULL; 316 const PetscInt *permutation_indices; 317 318 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 319 if (!is_simplex) { 320 PetscSection section; 321 322 // -- Get permutation 323 PetscCall(DMGetLocalSection(dm, §ion)); 324 PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation)); 325 PetscCall(ISGetIndices(permutation, &permutation_indices)); 326 } 327 328 // -- Copy interp, grad matrices 329 PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp)); 330 PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad)); 331 const CeedInt c = 0; 332 for (CeedInt q = 0; q < Q; q++) { 333 for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) { 334 CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp]; 335 336 interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c]; 337 for (CeedInt d = 0; d < dim; d++) { 338 grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d]; 339 } 340 } 341 } 342 343 // -- Cleanup 344 if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices)); 345 PetscCall(ISDestroy(&permutation)); 346 } 347 348 // Finally, create libCEED basis 349 PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis)); 350 PetscCall(PetscFree(q_points)); 351 PetscCall(PetscFree(interp)); 352 PetscCall(PetscFree(grad)); 353 PetscFunctionReturn(PETSC_SUCCESS); 354 } 355 356 // ----------------------------------------------------------------------------- 357 // Get CEED Basis from DMPlex 358 // ----------------------------------------------------------------------------- 359 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) { 360 PetscDS ds; 361 PetscFE fe; 362 PetscQuadrature quadrature; 363 PetscBool is_simplex = PETSC_TRUE; 364 PetscInt ds_field = -1; 365 366 PetscFunctionBeginUser; 367 // Get element information 368 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 369 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 370 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 371 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 372 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 373 374 // Check if simplex or tensor-product mesh 375 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 376 377 // Build libCEED basis 378 if (is_simplex) { 379 PetscTabulation basis_tabulation; 380 PetscInt num_derivatives = 1, face = 0; 381 382 PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation)); 383 PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis)); 384 } else { 385 PetscDualSpace dual_space; 386 PetscInt num_dual_basis_vectors; 387 PetscInt dim, num_comp, P, Q; 388 389 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 390 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 391 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 392 PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 393 P = num_dual_basis_vectors / num_comp; 394 PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL)); 395 396 CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim)); 397 CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim)); 398 399 PetscCallCeed(ceed, CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis)); 400 } 401 PetscFunctionReturn(PETSC_SUCCESS); 402 } 403 404 /** 405 @brief Setup `DM` with FE space of appropriate degree 406 407 Must be followed by `DMSetupByOrderEnd_FEM` 408 409 @param[in] setup_faces Flag to setup face geometry 410 @param[in] setup_coords Flag to setup coordinate spaces 411 @param[in] degree Polynomial orders of field 412 @param[in] coord_order Polynomial order of coordinate basis, or `RATEL_DECIDE` for default 413 @param[in] q_extra Additional quadrature order, or `RATEL_DECIDE` for default 414 @param[in] num_fields Number of fields in solution vector 415 @param[in] field_sizes Array of number of components for each field 416 @param[out] dm `DM` to setup 417 418 @return An error code: 0 - success, otherwise - failure 419 **/ 420 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 421 CeedInt num_fields, const CeedInt *field_sizes, DM dm) { 422 PetscInt dim, q_order = degree + q_extra; 423 PetscBool is_simplex = PETSC_TRUE; 424 PetscFE fe; 425 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 426 427 PetscFunctionBeginUser; 428 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 429 430 // Setup DM 431 PetscCall(DMGetDimension(dm, &dim)); 432 for (PetscInt i = 0; i < num_fields; i++) { 433 PetscFE fe_face; 434 PetscInt q_order = degree + q_extra; 435 436 PetscCall(PetscFECreateLagrange(comm, dim, field_sizes[i], is_simplex, degree, q_order, &fe)); 437 if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe_face)); 438 PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 439 PetscCall(PetscFEDestroy(&fe)); 440 } 441 PetscCall(DMCreateDS(dm)); 442 443 // Project coordinates to enrich quadrature space 444 if (setup_coords) { 445 DM dm_coord; 446 PetscDS ds_coord; 447 PetscFE fe_coord_current, fe_coord_new, fe_coord_face_new; 448 PetscDualSpace fe_coord_dual_space; 449 PetscInt fe_coord_order, num_comp_coord; 450 451 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 452 PetscCall(DMGetCoordinateDim(dm, &num_comp_coord)); 453 PetscCall(DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord, NULL)); 454 PetscCall(PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current)); 455 PetscCall(PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space)); 456 PetscCall(PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order)); 457 458 // Create FE for coordinates 459 PetscCheck(fe_coord_order == 1 || coord_order == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER_INPUT, 460 "Only linear mesh geometry supported. Recieved %d order", fe_coord_order); 461 PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, q_order, &fe_coord_new)); 462 if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe_coord_new, 1, &fe_coord_face_new)); 463 PetscCall(DMProjectCoordinates(dm, fe_coord_new)); 464 PetscCall(PetscFEDestroy(&fe_coord_new)); 465 } 466 PetscFunctionReturn(PETSC_SUCCESS); 467 } 468 469 /** 470 @brief Finish setting up `DM` 471 472 Must be called after `DMSetupByOrderBegin_FEM` and all strong BCs have be declared via `DMAddBoundaries`. 473 474 @param[in] setup_coords Flag to setup coordinate spaces 475 @param[out] dm `DM` to setup 476 477 @return An error code: 0 - success, otherwise - failure 478 **/ 479 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm) { 480 PetscBool is_simplex; 481 482 PetscFunctionBeginUser; 483 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 484 // Set tensor permutation if needed 485 if (!is_simplex) { 486 PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 487 if (setup_coords) { 488 DM dm_coord; 489 490 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 491 PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 492 } 493 } 494 PetscFunctionReturn(PETSC_SUCCESS); 495 } 496 497 /** 498 @brief Setup `DM` with FE space of appropriate degree with no boundary conditions 499 500 Calls `DMSetupByOrderBegin_FEM` and `DMSetupByOrderEnd_FEM` successively 501 502 @param[in] setup_faces Flag to setup face geometry 503 @param[in] setup_coords Flag to setup coordinate spaces 504 @param[in] degree Polynomial orders of field 505 @param[in] coord_order Polynomial order of coordinate basis, or `RATEL_DECIDE` for default 506 @param[in] q_extra Additional quadrature order, or `RATEL_DECIDE` for default 507 @param[in] num_fields Number of fields in solution vector 508 @param[in] field_sizes Array of number of components for each field 509 @param[out] dm `DM` to setup 510 511 @return An error code: 0 - success, otherwise - failure 512 **/ 513 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 514 CeedInt num_fields, const CeedInt *field_sizes, DM dm) { 515 PetscFunctionBeginUser; 516 PetscCall(DMSetupByOrderBegin_FEM(setup_faces, setup_coords, degree, coord_order, q_extra, num_fields, field_sizes, dm)); 517 PetscCall(DMSetupByOrderEnd_FEM(setup_coords, dm)); 518 PetscFunctionReturn(PETSC_SUCCESS); 519 } 520