17b8c5038SJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 27b8c5038SJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 37b8c5038SJames Wright 47b8c5038SJames Wright /// @file 57b8c5038SJames Wright /// Utilities for setting up DM and PetscFE 67b8c5038SJames Wright 77b8c5038SJames Wright #include <ceed.h> 87b8c5038SJames Wright #include <dm-utils.h> 97b8c5038SJames Wright #include <petsc-ceed-utils.h> 107b8c5038SJames Wright #include <petscdmplex.h> 117b8c5038SJames Wright #include <petscds.h> 127b8c5038SJames Wright 137b8c5038SJames Wright /** 147b8c5038SJames Wright @brief Convert `DM` field to `DS` field. 157b8c5038SJames Wright 167b8c5038SJames Wright @param[in] dm `DM` holding mesh 177b8c5038SJames Wright @param[in] domain_label Label for `DM` domain 187b8c5038SJames Wright @param[in] dm_field Index of `DM` field 197b8c5038SJames Wright @param[out] ds_field Index of `DS` field 207b8c5038SJames Wright 217b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 227b8c5038SJames Wright **/ 237b8c5038SJames Wright PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) { 247b8c5038SJames Wright PetscDS ds; 257b8c5038SJames Wright IS field_is; 267b8c5038SJames Wright const PetscInt *fields; 277b8c5038SJames Wright PetscInt num_fields; 287b8c5038SJames Wright 297b8c5038SJames Wright PetscFunctionBeginUser; 307b8c5038SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL)); 317b8c5038SJames Wright PetscCall(ISGetIndices(field_is, &fields)); 327b8c5038SJames Wright PetscCall(ISGetSize(field_is, &num_fields)); 337b8c5038SJames Wright for (PetscInt i = 0; i < num_fields; i++) { 347b8c5038SJames Wright if (dm_field == fields[i]) { 357b8c5038SJames Wright *ds_field = i; 367b8c5038SJames Wright break; 377b8c5038SJames Wright } 387b8c5038SJames Wright } 397b8c5038SJames Wright PetscCall(ISRestoreIndices(field_is, &fields)); 407b8c5038SJames Wright 417b8c5038SJames Wright PetscCheck(*ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 427b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 437b8c5038SJames Wright } 447b8c5038SJames Wright 457b8c5038SJames Wright /** 467b8c5038SJames Wright @brief Destroy `ElemRestriction` in a `PetscContainer`. 477b8c5038SJames Wright 487b8c5038SJames Wright Not collective across MPI processes. 497b8c5038SJames Wright 507b8c5038SJames Wright @param[in,out] ctx `CeedElemRestriction` 517b8c5038SJames Wright 527b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 537b8c5038SJames Wright **/ 547b8c5038SJames Wright static PetscErrorCode DMPlexCeedElemRestrictionDestroy(void **ctx) { 557b8c5038SJames Wright CeedElemRestriction rstr = *(CeedElemRestriction *)ctx; 567b8c5038SJames Wright 577b8c5038SJames Wright PetscFunctionBegin; 587b8c5038SJames Wright PetscCallCeed(CeedElemRestrictionReturnCeed(rstr), CeedElemRestrictionDestroy(&rstr)); 597b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 607b8c5038SJames Wright } 617b8c5038SJames Wright 627b8c5038SJames Wright /** 637b8c5038SJames Wright @brief Create `CeedElemRestriction` from `DMPlex`. 647b8c5038SJames Wright 657b8c5038SJames Wright The `CeedElemRestriction` is stored in the `DMPlex` object for reuse; 667b8c5038SJames Wright if the routine is called again with the same arguments, the same `CeedElemRestriction` is returned. 677b8c5038SJames Wright 687b8c5038SJames Wright Not collective across MPI processes. 697b8c5038SJames Wright 707b8c5038SJames Wright @param[in] ceed `Ceed` context 717b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 727b8c5038SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 737b8c5038SJames Wright @param[in] label_value Stratum value 747b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 757b8c5038SJames Wright @param[in] dm_field Index of `DMPlex` field 767b8c5038SJames Wright @param[out] restriction `CeedElemRestriction` for `DMPlex` 777b8c5038SJames Wright 787b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 797b8c5038SJames Wright **/ 807b8c5038SJames Wright PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 817b8c5038SJames Wright CeedElemRestriction *restriction) { 827b8c5038SJames Wright char container_name[1024]; 837b8c5038SJames Wright CeedElemRestriction container_restriction; 847b8c5038SJames Wright 857b8c5038SJames Wright PetscFunctionBeginUser; 867b8c5038SJames Wright { 877b8c5038SJames Wright const char *label_name = NULL; 887b8c5038SJames Wright 897b8c5038SJames Wright if (domain_label) PetscCall(PetscObjectGetName((PetscObject)domain_label, &label_name)); 907b8c5038SJames Wright PetscCall(PetscSNPrintf(container_name, sizeof(container_name), 917b8c5038SJames Wright "DM CeedElemRestriction - label: %s label value: %" PetscInt_FMT " height: %" PetscInt_FMT " DM field: %" PetscInt_FMT, 927b8c5038SJames Wright label_name ? label_name : "NONE", label_value, height, dm_field)); 937b8c5038SJames Wright } 947b8c5038SJames Wright PetscCall(PetscObjectContainerQuery((PetscObject)dm, container_name, (void **)&container_restriction)); 957b8c5038SJames Wright 967b8c5038SJames Wright if (container_restriction) { 977b8c5038SJames Wright // Copy existing restriction 987b8c5038SJames Wright *restriction = NULL; 997b8c5038SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(container_restriction, restriction)); 1007b8c5038SJames Wright } else { 1017b8c5038SJames Wright PetscInt num_elem, elem_size, num_dof, num_comp, *restriction_offsets_petsc; 1027b8c5038SJames Wright CeedInt *restriction_offsets_ceed = NULL; 1037b8c5038SJames Wright 1047b8c5038SJames Wright // Build restriction 1057b8c5038SJames Wright PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, 1067b8c5038SJames Wright &restriction_offsets_petsc)); 1077b8c5038SJames Wright PetscCall(IntArrayPetscToCeed(num_elem * elem_size, &restriction_offsets_petsc, &restriction_offsets_ceed)); 1087b8c5038SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 1097b8c5038SJames Wright restriction_offsets_ceed, restriction)); 1107b8c5038SJames Wright PetscCall(PetscFree(restriction_offsets_ceed)); 1117b8c5038SJames Wright 1127b8c5038SJames Wright // Set in container 1137b8c5038SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(*restriction, &container_restriction)); 1147b8c5038SJames Wright PetscCall(PetscObjectContainerCompose((PetscObject)dm, container_name, container_restriction, DMPlexCeedElemRestrictionDestroy)); 1157b8c5038SJames Wright } 1167b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1177b8c5038SJames Wright } 1187b8c5038SJames Wright 1197b8c5038SJames Wright /** 1207b8c5038SJames Wright @brief Create `CeedElemRestriction` from `DMPlex` domain for mesh coordinates. 1217b8c5038SJames Wright 1227b8c5038SJames Wright The `CeedElemRestriction` is stored in the coordinate `DMPlex` object for reuse; 1237b8c5038SJames Wright if the routine is called again with the same arguments, the same `CeedElemRestriction` is returned. 1247b8c5038SJames Wright 1257b8c5038SJames Wright Not collective across MPI processes. 1267b8c5038SJames Wright 1277b8c5038SJames Wright @param[in] ceed `Ceed` context 1287b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 1297b8c5038SJames Wright @param[in] domain_label Label for `DMPlex` domain 1307b8c5038SJames Wright @param[in] label_value Stratum value 1317b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 1327b8c5038SJames Wright @param[out] restriction `CeedElemRestriction` for mesh 1337b8c5038SJames Wright 1347b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 1357b8c5038SJames Wright **/ 1367b8c5038SJames Wright PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 1377b8c5038SJames Wright CeedElemRestriction *restriction) { 1387b8c5038SJames Wright DM dm_coord; 1397b8c5038SJames Wright 1407b8c5038SJames Wright PetscFunctionBeginUser; 1417b8c5038SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 1427b8c5038SJames Wright if (!dm_coord) { 1437b8c5038SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1447b8c5038SJames Wright } 1457b8c5038SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm_coord, domain_label, label_value, height, 0, restriction)); 1467b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1477b8c5038SJames Wright } 1487b8c5038SJames Wright 1497b8c5038SJames Wright /** 1507b8c5038SJames Wright @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data. 1517b8c5038SJames Wright 1527b8c5038SJames Wright Not collective across MPI processes. 1537b8c5038SJames Wright 1547b8c5038SJames Wright @param[in] ceed `Ceed` context 1557b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 1567b8c5038SJames Wright @param[in] domain_label Label for `DMPlex` domain 1577b8c5038SJames Wright @param[in] label_value Stratum value 1587b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 1597b8c5038SJames Wright @param[in] q_data_size Number of components for `QFunction` data 1607b8c5038SJames Wright @param[in] is_collocated Boolean flag indicating if the data is collocated on the nodes (`PETSC_TRUE`) or on quadrature points (`PETSC_FALSE`) 1617b8c5038SJames Wright @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 1627b8c5038SJames Wright 1637b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 1647b8c5038SJames Wright **/ 1657b8c5038SJames Wright static PetscErrorCode DMPlexCeedElemRestrictionStridedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 1667b8c5038SJames Wright PetscInt q_data_size, PetscBool is_collocated, CeedElemRestriction *restriction) { 1677b8c5038SJames Wright PetscInt num_elem, num_qpts, dm_field = 0; 1687b8c5038SJames Wright 1697b8c5038SJames Wright PetscFunctionBeginUser; 1707b8c5038SJames Wright { // Get number of elements 1717b8c5038SJames Wright PetscInt depth; 1727b8c5038SJames Wright DMLabel depth_label; 1737b8c5038SJames Wright IS point_is, depth_is; 1747b8c5038SJames Wright 1757b8c5038SJames Wright PetscCall(DMPlexGetDepth(dm, &depth)); 1767b8c5038SJames Wright PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 1777b8c5038SJames Wright PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 1787b8c5038SJames Wright if (domain_label) { 1797b8c5038SJames Wright IS domain_is; 1807b8c5038SJames Wright 1817b8c5038SJames Wright PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is)); 1827b8c5038SJames Wright if (domain_is) { 1837b8c5038SJames Wright PetscCall(ISIntersect(depth_is, domain_is, &point_is)); 1847b8c5038SJames Wright PetscCall(ISDestroy(&domain_is)); 1857b8c5038SJames Wright } else { 1867b8c5038SJames Wright point_is = NULL; 1877b8c5038SJames Wright } 1887b8c5038SJames Wright PetscCall(ISDestroy(&depth_is)); 1897b8c5038SJames Wright } else { 1907b8c5038SJames Wright point_is = depth_is; 1917b8c5038SJames Wright } 1927b8c5038SJames Wright if (point_is) { 1937b8c5038SJames Wright PetscCall(ISGetLocalSize(point_is, &num_elem)); 1947b8c5038SJames Wright } else { 1957b8c5038SJames Wright num_elem = 0; 1967b8c5038SJames Wright } 1977b8c5038SJames Wright PetscCall(ISDestroy(&point_is)); 1987b8c5038SJames Wright } 1997b8c5038SJames Wright 2007b8c5038SJames Wright { // Get number of quadrature points 2017b8c5038SJames Wright PetscDS ds; 2027b8c5038SJames Wright PetscFE fe; 2037b8c5038SJames Wright PetscInt ds_field = -1; 2047b8c5038SJames Wright 2057b8c5038SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 2067b8c5038SJames Wright PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 2077b8c5038SJames Wright PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 2087b8c5038SJames Wright PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 2097b8c5038SJames Wright if (is_collocated) { 2107b8c5038SJames Wright PetscDualSpace dual_space; 2117b8c5038SJames Wright PetscInt num_dual_basis_vectors, dim, num_comp; 2127b8c5038SJames Wright 2137b8c5038SJames Wright PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 2147b8c5038SJames Wright PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 2157b8c5038SJames Wright PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 2167b8c5038SJames Wright PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 2177b8c5038SJames Wright num_qpts = num_dual_basis_vectors / num_comp; 2187b8c5038SJames Wright } else { 2197b8c5038SJames Wright PetscQuadrature quadrature; 2207b8c5038SJames Wright 2217b8c5038SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 2227b8c5038SJames Wright PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 2237b8c5038SJames Wright PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 2247b8c5038SJames Wright PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 2257b8c5038SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 2267b8c5038SJames Wright PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &num_qpts, NULL, NULL)); 2277b8c5038SJames Wright } 2287b8c5038SJames Wright } 2297b8c5038SJames Wright 2307b8c5038SJames Wright // Create the restriction 2317b8c5038SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, num_elem * num_qpts * q_data_size, CEED_STRIDES_BACKEND, 2327b8c5038SJames Wright restriction)); 2337b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2347b8c5038SJames Wright } 2357b8c5038SJames Wright 2367b8c5038SJames Wright /** 2377b8c5038SJames Wright @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data. 2387b8c5038SJames Wright 2397b8c5038SJames Wright Not collective across MPI processes. 2407b8c5038SJames Wright 2417b8c5038SJames Wright @param[in] ceed `Ceed` context 2427b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 2437b8c5038SJames Wright @param[in] domain_label Label for `DMPlex` domain 2447b8c5038SJames Wright @param[in] label_value Stratum value 2457b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 2467b8c5038SJames Wright @param[in] q_data_size Number of components for `QFunction` data 2477b8c5038SJames Wright @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 2487b8c5038SJames Wright 2497b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 2507b8c5038SJames Wright **/ 2517b8c5038SJames Wright PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 2527b8c5038SJames Wright PetscInt q_data_size, CeedElemRestriction *restriction) { 2537b8c5038SJames Wright PetscFunctionBeginUser; 2547b8c5038SJames Wright PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_FALSE, restriction)); 2557b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2567b8c5038SJames Wright } 2577b8c5038SJames Wright 2587b8c5038SJames Wright /** 2597b8c5038SJames Wright @brief Create `CeedElemRestriction` from `DMPlex` domain for nodally collocated auxilury `QFunction` data. 2607b8c5038SJames Wright 2617b8c5038SJames Wright Not collective across MPI processes. 2627b8c5038SJames Wright 2637b8c5038SJames Wright @param[in] ceed `Ceed` context 2647b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 2657b8c5038SJames Wright @param[in] domain_label Label for `DMPlex` domain 2667b8c5038SJames Wright @param[in] label_value Stratum value 2677b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 2687b8c5038SJames Wright @param[in] q_data_size Number of components for `QFunction` data 2697b8c5038SJames Wright @param[out] restriction Strided `CeedElemRestriction` for `QFunction` data 2707b8c5038SJames Wright 2717b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 2727b8c5038SJames Wright **/ 2737b8c5038SJames Wright PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 2747b8c5038SJames Wright PetscInt q_data_size, CeedElemRestriction *restriction) { 2757b8c5038SJames Wright PetscFunctionBeginUser; 2767b8c5038SJames Wright PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_TRUE, restriction)); 2777b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2787b8c5038SJames Wright } 2797b8c5038SJames Wright 2807b8c5038SJames Wright /** 2817b8c5038SJames Wright @brief Creates a tensor-product uniform quadrature. 2827b8c5038SJames Wright This is only for comparison studies and generally should not be used. 2837b8c5038SJames Wright 2847b8c5038SJames Wright @param[in] dim Spatial dimension 2857b8c5038SJames Wright @param[in] num_comp Number of components 2867b8c5038SJames Wright @param[in] num_points Number of points in one dimension 2877b8c5038SJames Wright @param[in] a Left end of interval (often -1) 2887b8c5038SJames Wright @param[in] b Right end of interval (often +1) 2897b8c5038SJames Wright @param[out] q `PetscQuadrature` object 2907b8c5038SJames Wright 2917b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 2927b8c5038SJames Wright **/ 2937b8c5038SJames Wright static PetscErrorCode PetscDTUniformTensorQuadrature(PetscInt dim, PetscInt num_comp, PetscInt num_points, PetscReal a, PetscReal b, 2947b8c5038SJames Wright PetscQuadrature *q) { 2957b8c5038SJames Wright DMPolytopeType ct; 2967b8c5038SJames Wright PetscInt num_total_points = dim > 1 ? (dim > 2 ? (num_points * num_points * num_points) : (num_points * num_points)) : num_points; 2977b8c5038SJames Wright PetscReal *coords, *weights, *coords_1d; 2987b8c5038SJames Wright 2997b8c5038SJames Wright PetscFunctionBeginUser; 3007b8c5038SJames Wright PetscCall(PetscMalloc1(num_total_points * dim, &coords)); 3017b8c5038SJames Wright PetscCall(PetscMalloc1(num_total_points * num_comp, &weights)); 3027b8c5038SJames Wright // Compute weights, points 3037b8c5038SJames Wright switch (dim) { 3047b8c5038SJames Wright case 0: 3057b8c5038SJames Wright ct = DM_POLYTOPE_POINT; 3067b8c5038SJames Wright PetscCall(PetscFree(coords)); 3077b8c5038SJames Wright PetscCall(PetscFree(weights)); 3087b8c5038SJames Wright PetscCall(PetscMalloc1(1, &coords)); 3097b8c5038SJames Wright PetscCall(PetscMalloc1(num_comp, &weights)); 3107b8c5038SJames Wright num_total_points = 1; 3117b8c5038SJames Wright coords[0] = 0.0; 3127b8c5038SJames Wright for (PetscInt c = 0; c < num_comp; c++) weights[c] = 1.0; 3137b8c5038SJames Wright break; 3147b8c5038SJames Wright case 1: { 3157b8c5038SJames Wright ct = DM_POLYTOPE_SEGMENT; 3167b8c5038SJames Wright PetscReal step = (b - a) / num_points; 3177b8c5038SJames Wright 3187b8c5038SJames Wright for (PetscInt i = 0; i < num_points; i++) { 3197b8c5038SJames Wright coords[i] = step * (i + 0.5) + a; 3207b8c5038SJames Wright for (PetscInt c = 0; c < num_comp; c++) weights[i * num_comp + c] = 1.0; 3217b8c5038SJames Wright } 3227b8c5038SJames Wright } break; 3237b8c5038SJames Wright case 2: { 3247b8c5038SJames Wright ct = DM_POLYTOPE_QUADRILATERAL; 3257b8c5038SJames Wright PetscCall(PetscMalloc1(num_points, &coords_1d)); 3267b8c5038SJames Wright PetscReal step = (b - a) / num_points; 3277b8c5038SJames Wright 3287b8c5038SJames Wright for (PetscInt i = 0; i < num_points; i++) coords_1d[i] = step * (i + 0.5) + a; 3297b8c5038SJames Wright for (PetscInt i = 0; i < num_points; i++) { 3307b8c5038SJames Wright for (PetscInt j = 0; j < num_points; j++) { 3317b8c5038SJames Wright coords[(i * num_points + j) * dim + 0] = coords_1d[i]; 3327b8c5038SJames Wright coords[(i * num_points + j) * dim + 1] = coords_1d[j]; 3337b8c5038SJames Wright for (PetscInt c = 0; c < num_comp; c++) weights[(i * num_points + j) * num_comp + c] = 1.0; 3347b8c5038SJames Wright } 3357b8c5038SJames Wright } 3367b8c5038SJames Wright PetscCall(PetscFree(coords_1d)); 3377b8c5038SJames Wright } break; 3387b8c5038SJames Wright case 3: { 3397b8c5038SJames Wright ct = DM_POLYTOPE_HEXAHEDRON; 3407b8c5038SJames Wright PetscCall(PetscMalloc1(num_points, &coords_1d)); 3417b8c5038SJames Wright PetscReal step = (b - a) / num_points; 3427b8c5038SJames Wright 3437b8c5038SJames Wright for (PetscInt i = 0; i < num_points; i++) coords_1d[i] = step * (i + 0.5) + a; 3447b8c5038SJames Wright for (PetscInt i = 0; i < num_points; i++) { 3457b8c5038SJames Wright for (PetscInt j = 0; j < num_points; j++) { 3467b8c5038SJames Wright for (PetscInt k = 0; k < num_points; k++) { 3477b8c5038SJames Wright coords[((i * num_points + j) * num_points + k) * dim + 0] = coords_1d[i]; 3487b8c5038SJames Wright coords[((i * num_points + j) * num_points + k) * dim + 1] = coords_1d[j]; 3497b8c5038SJames Wright coords[((i * num_points + j) * num_points + k) * dim + 2] = coords_1d[k]; 3507b8c5038SJames Wright for (PetscInt c = 0; c < num_comp; c++) weights[((i * num_points + j) * num_points + k) * num_comp + c] = 1.0; 3517b8c5038SJames Wright } 3527b8c5038SJames Wright } 3537b8c5038SJames Wright } 3547b8c5038SJames Wright PetscCall(PetscFree(coords_1d)); 3557b8c5038SJames Wright } break; 3567b8c5038SJames Wright // LCOV_EXCL_START 3577b8c5038SJames Wright default: 3587b8c5038SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %" PetscInt_FMT, dim); 3597b8c5038SJames Wright // LCOV_EXCL_STOP 3607b8c5038SJames Wright } 3617b8c5038SJames Wright PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, q)); 3627b8c5038SJames Wright PetscCall(PetscQuadratureSetCellType(*q, ct)); 3637b8c5038SJames Wright PetscCall(PetscQuadratureSetOrder(*q, num_points)); 3647b8c5038SJames Wright PetscCall(PetscQuadratureSetData(*q, dim, num_comp, num_total_points, coords, weights)); 3657b8c5038SJames Wright PetscCall(PetscObjectChangeTypeName((PetscObject)*q, "UniformTensor")); 3667b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3677b8c5038SJames Wright } 3687b8c5038SJames Wright 3697b8c5038SJames Wright /** 3707b8c5038SJames Wright @brief Setup `DM` with FE space of appropriate degree 3717b8c5038SJames Wright 3727b8c5038SJames Wright @param[in] comm MPI communicator 3737b8c5038SJames Wright @param[in] dim Spatial dimension 3747b8c5038SJames Wright @param[in] num_comp Number of components 3757b8c5038SJames Wright @param[in] is_simplex Flag for simplex or tensor product element 3767b8c5038SJames Wright @param[in] order Polynomial order of space 3777b8c5038SJames Wright @param[in] q_order Quadrature order 3787b8c5038SJames Wright @param[in] prefix The options prefix, or `NULL` 3797b8c5038SJames Wright @param[out] fem `PetscFE` object 3807b8c5038SJames Wright 3817b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 3827b8c5038SJames Wright **/ 3837b8c5038SJames Wright PetscErrorCode PetscFECreateLagrangeFromOptions(MPI_Comm comm, PetscInt dim, PetscInt num_comp, PetscBool is_simplex, PetscInt order, 3847b8c5038SJames Wright PetscInt q_order, const char prefix[], PetscFE *fem) { 3857b8c5038SJames Wright PetscBool is_tensor = !is_simplex; 3867b8c5038SJames Wright DMPolytopeType polytope_type = DMPolytopeTypeSimpleShape(dim, is_simplex); 3877b8c5038SJames Wright PetscSpace fe_space; 3887b8c5038SJames Wright PetscDualSpace fe_dual_space; 3897b8c5038SJames Wright PetscQuadrature quadrature, face_quadrature; 3907b8c5038SJames Wright 3917b8c5038SJames Wright PetscFunctionBeginUser; 3927b8c5038SJames Wright // Create space 3937b8c5038SJames Wright PetscCall(PetscSpaceCreate(comm, &fe_space)); 3947b8c5038SJames Wright PetscCall(PetscSpaceSetType(fe_space, PETSCSPACEPOLYNOMIAL)); 3957b8c5038SJames Wright PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fe_space, prefix)); 3967b8c5038SJames Wright PetscCall(PetscSpacePolynomialSetTensor(fe_space, is_tensor)); 3977b8c5038SJames Wright PetscCall(PetscSpaceSetNumComponents(fe_space, num_comp)); 3987b8c5038SJames Wright PetscCall(PetscSpaceSetNumVariables(fe_space, dim)); 3997b8c5038SJames Wright if (order >= 0) { 4007b8c5038SJames Wright PetscCall(PetscSpaceSetDegree(fe_space, order, PETSC_DETERMINE)); 4017b8c5038SJames Wright if (polytope_type == DM_POLYTOPE_TRI_PRISM || polytope_type == DM_POLYTOPE_TRI_PRISM_TENSOR) { 4027b8c5038SJames Wright PetscSpace fe_space_end, fe_space_side; 4037b8c5038SJames Wright 4047b8c5038SJames Wright PetscCall(PetscSpaceSetNumComponents(fe_space, 1)); 4057b8c5038SJames Wright PetscCall(PetscSpaceCreate(comm, &fe_space_end)); 4067b8c5038SJames Wright PetscCall(PetscSpaceSetType(fe_space_end, PETSCSPACEPOLYNOMIAL)); 4077b8c5038SJames Wright PetscCall(PetscSpacePolynomialSetTensor(fe_space_end, PETSC_FALSE)); 4087b8c5038SJames Wright PetscCall(PetscSpaceSetNumComponents(fe_space_end, 1)); 4097b8c5038SJames Wright PetscCall(PetscSpaceSetNumVariables(fe_space_end, dim - 1)); 4107b8c5038SJames Wright PetscCall(PetscSpaceSetDegree(fe_space_end, order, PETSC_DETERMINE)); 4117b8c5038SJames Wright PetscCall(PetscSpaceCreate(comm, &fe_space_side)); 4127b8c5038SJames Wright PetscCall(PetscSpaceSetType(fe_space_side, PETSCSPACEPOLYNOMIAL)); 4137b8c5038SJames Wright PetscCall(PetscSpacePolynomialSetTensor(fe_space_side, PETSC_FALSE)); 4147b8c5038SJames Wright PetscCall(PetscSpaceSetNumComponents(fe_space_side, 1)); 4157b8c5038SJames Wright PetscCall(PetscSpaceSetNumVariables(fe_space_side, 1)); 4167b8c5038SJames Wright PetscCall(PetscSpaceSetDegree(fe_space_side, order, PETSC_DETERMINE)); 4177b8c5038SJames Wright PetscCall(PetscSpaceSetType(fe_space, PETSCSPACETENSOR)); 4187b8c5038SJames Wright PetscCall(PetscSpaceTensorSetNumSubspaces(fe_space, 2)); 4197b8c5038SJames Wright PetscCall(PetscSpaceTensorSetSubspace(fe_space, 0, fe_space_end)); 4207b8c5038SJames Wright PetscCall(PetscSpaceTensorSetSubspace(fe_space, 1, fe_space_side)); 4217b8c5038SJames Wright PetscCall(PetscSpaceDestroy(&fe_space_end)); 4227b8c5038SJames Wright PetscCall(PetscSpaceDestroy(&fe_space_side)); 4237b8c5038SJames Wright 4247b8c5038SJames Wright if (num_comp > 1) { 4257b8c5038SJames Wright PetscSpace scalar_fe_space = fe_space; 4267b8c5038SJames Wright 4277b8c5038SJames Wright PetscCall(PetscSpaceCreate(comm, &fe_space)); 4287b8c5038SJames Wright PetscCall(PetscSpaceSetNumVariables(fe_space, dim)); 4297b8c5038SJames Wright PetscCall(PetscSpaceSetNumComponents(fe_space, num_comp)); 4307b8c5038SJames Wright PetscCall(PetscSpaceSetType(fe_space, PETSCSPACESUM)); 4317b8c5038SJames Wright PetscCall(PetscSpaceSumSetNumSubspaces(fe_space, num_comp)); 4327b8c5038SJames Wright PetscCall(PetscSpaceSumSetConcatenate(fe_space, PETSC_TRUE)); 4337b8c5038SJames Wright PetscCall(PetscSpaceSumSetInterleave(fe_space, PETSC_TRUE, PETSC_FALSE)); 4347b8c5038SJames Wright for (PetscInt i = 0; i < num_comp; i++) PetscCall(PetscSpaceSumSetSubspace(fe_space, i, scalar_fe_space)); 4357b8c5038SJames Wright PetscCall(PetscSpaceDestroy(&scalar_fe_space)); 4367b8c5038SJames Wright } 4377b8c5038SJames Wright } 4387b8c5038SJames Wright } 4397b8c5038SJames Wright PetscCall(PetscSpaceSetFromOptions(fe_space)); 4407b8c5038SJames Wright PetscCall(PetscSpaceSetUp(fe_space)); 4417b8c5038SJames Wright PetscCall(PetscSpaceGetDegree(fe_space, &order, NULL)); 4427b8c5038SJames Wright PetscCall(PetscSpacePolynomialGetTensor(fe_space, &is_tensor)); 4437b8c5038SJames Wright PetscCall(PetscSpaceGetNumComponents(fe_space, &num_comp)); 4447b8c5038SJames Wright 4457b8c5038SJames Wright // Create dual space 4467b8c5038SJames Wright PetscCall(PetscDualSpaceCreate(comm, &fe_dual_space)); 4477b8c5038SJames Wright PetscCall(PetscDualSpaceSetType(fe_dual_space, PETSCDUALSPACELAGRANGE)); 4487b8c5038SJames Wright PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fe_dual_space, prefix)); 4497b8c5038SJames Wright { 4507b8c5038SJames Wright DM dual_space_dm; 4517b8c5038SJames Wright 4527b8c5038SJames Wright PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, polytope_type, &dual_space_dm)); 4537b8c5038SJames Wright PetscCall(PetscDualSpaceSetDM(fe_dual_space, dual_space_dm)); 4547b8c5038SJames Wright PetscCall(DMDestroy(&dual_space_dm)); 4557b8c5038SJames Wright } 4567b8c5038SJames Wright PetscCall(PetscDualSpaceSetNumComponents(fe_dual_space, num_comp)); 4577b8c5038SJames Wright PetscCall(PetscDualSpaceSetOrder(fe_dual_space, order)); 4587b8c5038SJames Wright PetscCall(PetscDualSpaceLagrangeSetTensor(fe_dual_space, (is_tensor || (polytope_type == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE)); 4597b8c5038SJames Wright PetscCall(PetscDualSpaceSetFromOptions(fe_dual_space)); 4607b8c5038SJames Wright PetscCall(PetscDualSpaceSetUp(fe_dual_space)); 4617b8c5038SJames Wright 4627b8c5038SJames Wright // Create quadrature 4637b8c5038SJames Wright q_order = q_order >= 0 ? q_order : order; 4647b8c5038SJames Wright PetscBool use_uniform = PETSC_FALSE; 4657b8c5038SJames Wright 4667b8c5038SJames Wright PetscOptionsBegin(comm, NULL, "Uniform quadrature check", NULL); 4677b8c5038SJames Wright PetscCall(PetscOptionsBool("-petscdt_uniform_tensor_quadrature", "", NULL, use_uniform, &use_uniform, NULL)); 4687b8c5038SJames Wright PetscOptionsEnd(); 4697b8c5038SJames Wright PetscCheck(!use_uniform || (is_tensor || dim == 1), comm, PETSC_ERR_SUP, "Can only use uniform quadrature on tensor product elements"); 4707b8c5038SJames Wright if (use_uniform) { 4717b8c5038SJames Wright PetscCall(PetscDTUniformTensorQuadrature(dim, 1, q_order + 1, -1.0, 1.0, &quadrature)); 4727b8c5038SJames Wright PetscCall(PetscDTUniformTensorQuadrature(dim - 1, 1, q_order + 1, -1.0, 1.0, &face_quadrature)); 4737b8c5038SJames Wright } else { 4747b8c5038SJames Wright PetscCall(PetscDTCreateDefaultQuadrature(polytope_type, q_order, &quadrature, &face_quadrature)); 4757b8c5038SJames Wright } 4767b8c5038SJames Wright 4777b8c5038SJames Wright // Create finite element 4787b8c5038SJames Wright PetscCall(PetscFECreateFromSpaces(fe_space, fe_dual_space, quadrature, face_quadrature, fem)); 4797b8c5038SJames Wright PetscCall(PetscFESetFromOptions(*fem)); 4807b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4817b8c5038SJames Wright } 4827b8c5038SJames Wright 4837b8c5038SJames Wright /** 4847b8c5038SJames Wright @brief Get global `DMPlex` topology type. 4857b8c5038SJames Wright 4867b8c5038SJames Wright Collective across MPI processes. 4877b8c5038SJames Wright 4887b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 4897b8c5038SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 4907b8c5038SJames Wright @param[in] label_value Stratum value 4917b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 4927b8c5038SJames Wright @param[out] cell_type `DMPlex` topology type 4937b8c5038SJames Wright **/ 4947b8c5038SJames Wright static inline PetscErrorCode GetGlobalDMPlexPolytopeType(DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 4957b8c5038SJames Wright DMPolytopeType *cell_type) { 4967b8c5038SJames Wright PetscInt first_point; 4977b8c5038SJames Wright PetscInt ids[1] = {label_value}; 4987b8c5038SJames Wright DMLabel depth_label; 4997b8c5038SJames Wright 5007b8c5038SJames Wright PetscFunctionBeginUser; 5017b8c5038SJames Wright // Use depth label if no domain label present 5027b8c5038SJames Wright if (!domain_label) { 5037b8c5038SJames Wright PetscInt depth; 5047b8c5038SJames Wright 5057b8c5038SJames Wright PetscCall(DMPlexGetDepth(dm, &depth)); 5067b8c5038SJames Wright PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 5077b8c5038SJames Wright ids[0] = depth - height; 5087b8c5038SJames Wright } 5097b8c5038SJames Wright 5107b8c5038SJames Wright // Get cell interp, grad, and quadrature data 5117b8c5038SJames Wright PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL)); 5127b8c5038SJames Wright if (first_point != -1) PetscCall(DMPlexGetCellType(dm, first_point, cell_type)); 5137b8c5038SJames Wright 5147b8c5038SJames Wright { // Form agreement about CellType 5157b8c5038SJames Wright PetscInt cell_type_local = -1, cell_type_global = -1; 5167b8c5038SJames Wright 5177b8c5038SJames Wright if (first_point != -1) cell_type_local = (PetscInt)*cell_type; 5187b8c5038SJames Wright PetscCallMPI(MPIU_Allreduce(&cell_type_local, &cell_type_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 5197b8c5038SJames Wright *cell_type = (DMPolytopeType)cell_type_global; 5207b8c5038SJames Wright } 5217b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5227b8c5038SJames Wright } 5237b8c5038SJames Wright 5247b8c5038SJames Wright /** 5257b8c5038SJames Wright @brief Get the permutation and field offset for a given depth. 5267b8c5038SJames Wright 5277b8c5038SJames Wright Not collective across MPI processes. 5287b8c5038SJames Wright 5297b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 5307b8c5038SJames Wright @param[in] depth Depth of `DMPlex` topology 5317b8c5038SJames Wright @param[in] field Index of `DMPlex` field 5327b8c5038SJames Wright @param[out] permutation Permutation for `DMPlex` field 5337b8c5038SJames Wright @param[out] field_offset Offset to `DMPlex` field 5347b8c5038SJames Wright **/ 5357b8c5038SJames Wright static inline PetscErrorCode GetClosurePermutationAndFieldOffsetAtDepth(DM dm, PetscInt depth, PetscInt field, IS *permutation, 5367b8c5038SJames Wright PetscInt *field_offset) { 5377b8c5038SJames Wright PetscBool is_field_continuous = PETSC_TRUE; 5387b8c5038SJames Wright PetscInt dim, num_fields, offset = 0, size = 0; 5397b8c5038SJames Wright PetscSection section; 5407b8c5038SJames Wright 5417b8c5038SJames Wright PetscFunctionBeginUser; 5427b8c5038SJames Wright *field_offset = 0; 5437b8c5038SJames Wright 5447b8c5038SJames Wright PetscCall(DMGetDimension(dm, &dim)); 5457b8c5038SJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 5467b8c5038SJames Wright PetscCall(PetscSectionGetNumFields(section, &num_fields)); 5477b8c5038SJames Wright // ---- Permutation size and offset to current field 5487b8c5038SJames Wright for (PetscInt f = 0; f < num_fields; f++) { 5497b8c5038SJames Wright PetscBool is_continuous; 5507b8c5038SJames Wright PetscInt num_components, num_dof_1d, dual_space_size, field_size; 5517b8c5038SJames Wright PetscObject obj; 5527b8c5038SJames Wright PetscFE fe = NULL; 5537b8c5038SJames Wright PetscDualSpace dual_space; 5547b8c5038SJames Wright 5557b8c5038SJames Wright PetscCall(PetscSectionGetFieldComponents(section, f, &num_components)); 5567b8c5038SJames Wright PetscCall(DMGetField(dm, f, NULL, &obj)); 5577b8c5038SJames Wright fe = (PetscFE)obj; 5587b8c5038SJames Wright PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 5597b8c5038SJames Wright PetscCall(PetscDualSpaceLagrangeGetContinuity(dual_space, &is_continuous)); 5607b8c5038SJames Wright if (f == field) is_field_continuous = is_continuous; 5617b8c5038SJames Wright if (!is_continuous) continue; 5627b8c5038SJames Wright PetscCall(PetscDualSpaceGetDimension(dual_space, &dual_space_size)); 5637b8c5038SJames Wright num_dof_1d = (PetscInt)PetscCeilReal(PetscPowReal(dual_space_size / num_components, 1.0 / dim)); 5647b8c5038SJames Wright field_size = PetscPowInt(num_dof_1d, depth) * num_components; 5657b8c5038SJames Wright if (f < field) offset += field_size; 5667b8c5038SJames Wright size += field_size; 5677b8c5038SJames Wright } 5687b8c5038SJames Wright 5697b8c5038SJames Wright if (is_field_continuous) { 5707b8c5038SJames Wright *field_offset = offset; 5717b8c5038SJames Wright PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, depth, size, permutation)); 5727b8c5038SJames Wright } 5737b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5747b8c5038SJames Wright } 5757b8c5038SJames Wright 5767b8c5038SJames Wright /** 5777b8c5038SJames Wright @brief Compute field tabulation from `PetscTabulation`. 5787b8c5038SJames Wright 5797b8c5038SJames Wright Not collective across MPI processes. 5807b8c5038SJames Wright 5817b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 5827b8c5038SJames Wright @param[in] field Index of `DMPlex` field 5837b8c5038SJames Wright @param[in] face Index of basis face, or 0 5847b8c5038SJames Wright @param[in] tabulation PETSc basis tabulation 5857b8c5038SJames Wright @param[out] interp Interpolation matrix in libCEED orientation 5867b8c5038SJames Wright @param[out] grad Gradient matrix in libCEED orientation 5877b8c5038SJames Wright 5887b8c5038SJames Wright @note `interp` and `grad` are allocated using `PetscCalloc1` and must be freed by the user. 5897b8c5038SJames Wright **/ 5907b8c5038SJames Wright static inline PetscErrorCode ComputeFieldTabulationP2C(DM dm, PetscInt field, PetscInt face, PetscTabulation tabulation, CeedScalar **interp, 5917b8c5038SJames Wright CeedScalar **grad) { 5927b8c5038SJames Wright PetscBool is_simplex, has_permutation = PETSC_FALSE; 5937b8c5038SJames Wright PetscInt field_offset = 0, num_comp, P, Q, dim; 5947b8c5038SJames Wright const PetscInt *permutation_indices; 5957b8c5038SJames Wright IS permutation = NULL; 5967b8c5038SJames Wright 5977b8c5038SJames Wright PetscFunctionBeginUser; 5987b8c5038SJames Wright num_comp = tabulation->Nc; 5997b8c5038SJames Wright P = tabulation->Nb / num_comp; 6007b8c5038SJames Wright Q = tabulation->Np; 6017b8c5038SJames Wright dim = tabulation->cdim; 6027b8c5038SJames Wright 6037b8c5038SJames Wright PetscCall(PetscCalloc1(P * Q, interp)); 6047b8c5038SJames Wright PetscCall(PetscCalloc1(P * Q * dim, grad)); 6057b8c5038SJames Wright 6067b8c5038SJames Wright PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 6077b8c5038SJames Wright if (!is_simplex) { 6087b8c5038SJames Wright PetscCall(GetClosurePermutationAndFieldOffsetAtDepth(dm, dim, field, &permutation, &field_offset)); 6097b8c5038SJames Wright has_permutation = !!permutation; 6107b8c5038SJames Wright if (has_permutation) PetscCall(ISGetIndices(permutation, &permutation_indices)); 6117b8c5038SJames Wright } 6127b8c5038SJames Wright 6137b8c5038SJames Wright const CeedInt c = 0; 6147b8c5038SJames Wright for (CeedInt q = 0; q < Q; q++) { 6157b8c5038SJames Wright for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) { 6167b8c5038SJames Wright CeedInt p_petsc = has_permutation ? (permutation_indices[field_offset + p_ceed * num_comp] - field_offset) : (p_ceed * num_comp); 6177b8c5038SJames Wright (*interp)[q * P + p_ceed] = tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c]; 6187b8c5038SJames Wright for (CeedInt d = 0; d < dim; d++) { 6197b8c5038SJames Wright (*grad)[(d * Q + q) * P + p_ceed] = tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d]; 6207b8c5038SJames Wright } 6217b8c5038SJames Wright } 6227b8c5038SJames Wright } 6237b8c5038SJames Wright 6247b8c5038SJames Wright if (has_permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices)); 6257b8c5038SJames Wright PetscCall(ISDestroy(&permutation)); 6267b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6277b8c5038SJames Wright } 6287b8c5038SJames Wright 6297b8c5038SJames Wright /** 6307b8c5038SJames Wright @brief Get quadrature data from `PetscQuadrature` for use with libCEED. 6317b8c5038SJames Wright 6327b8c5038SJames Wright Not collective across MPI processes. 6337b8c5038SJames Wright 6347b8c5038SJames Wright @param[in] fe `PetscFE` object 6357b8c5038SJames Wright @param[in] quadrature PETSc basis quadrature 6367b8c5038SJames Wright @param[out] q_dim Quadrature dimension, or NULL if not needed 6377b8c5038SJames Wright @param[out] num_comp Number of components, or NULL if not needed 6387b8c5038SJames Wright @param[out] Q Number of quadrature points, or NULL if not needed 6397b8c5038SJames Wright @param[out] q_points Quadrature points in libCEED orientation, or NULL if not needed 6407b8c5038SJames Wright @param[out] q_weights Quadrature weights in libCEED orientation, or NULL if not needed 6417b8c5038SJames Wright 6427b8c5038SJames Wright @note `q_weights` and `q_points` are allocated using `PetscCalloc1` and must be freed by the user. 6437b8c5038SJames Wright **/ 6447b8c5038SJames Wright static inline PetscErrorCode GetQuadratureDataP2C(PetscFE fe, PetscQuadrature quadrature, PetscInt *q_dim, PetscInt *num_comp, PetscInt *Q, 6457b8c5038SJames Wright CeedScalar **q_points, CeedScalar **q_weights) { 6467b8c5038SJames Wright PetscInt spatial_dim, dim, num_comp_quadrature, num_q_points; 6477b8c5038SJames Wright const PetscReal *q_points_petsc, *q_weights_petsc; 6487b8c5038SJames Wright 6497b8c5038SJames Wright PetscFunctionBeginUser; 6507b8c5038SJames Wright PetscCall(PetscFEGetSpatialDimension(fe, &spatial_dim)); 6517b8c5038SJames Wright PetscCall(PetscQuadratureGetData(quadrature, &dim, &num_comp_quadrature, &num_q_points, &q_points_petsc, &q_weights_petsc)); 6527b8c5038SJames Wright if (q_dim) *q_dim = dim; 6537b8c5038SJames Wright if (num_comp) *num_comp = num_comp_quadrature; 6547b8c5038SJames Wright if (Q) *Q = num_q_points; 6557b8c5038SJames Wright if (q_weights) { 6567b8c5038SJames Wright PetscSizeT q_weights_size = num_q_points; 6577b8c5038SJames Wright 6587b8c5038SJames Wright PetscCall(PetscCalloc1(q_weights_size, q_weights)); 6597b8c5038SJames Wright for (CeedInt i = 0; i < num_q_points; i++) (*q_weights)[i] = q_weights_petsc[i]; 6607b8c5038SJames Wright } 6617b8c5038SJames Wright if (q_points) { 6627b8c5038SJames Wright PetscSizeT q_points_size = num_q_points * spatial_dim; 6637b8c5038SJames Wright 6647b8c5038SJames Wright PetscCall(PetscCalloc1(q_points_size, q_points)); 6657b8c5038SJames Wright for (CeedInt i = 0; i < num_q_points; i++) { 6667b8c5038SJames Wright for (CeedInt d = 0; d < dim; d++) (*q_points)[i * spatial_dim + d] = q_points_petsc[i * dim + d]; 6677b8c5038SJames Wright } 6687b8c5038SJames Wright } 6697b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6707b8c5038SJames Wright } 6717b8c5038SJames Wright 6727b8c5038SJames Wright /** 6737b8c5038SJames Wright @brief Create 1D tabulation from `PetscFE`. 6747b8c5038SJames Wright 6757b8c5038SJames Wright Collective across MPI processes. 6767b8c5038SJames Wright 6777b8c5038SJames Wright @note `q_weights` and `q_points` are allocated using `PetscCalloc1` and must be freed by the user. 6787b8c5038SJames Wright 6797b8c5038SJames Wright @param[in] comm `MPI_Comm` for creating `FE` object from options 6807b8c5038SJames Wright @param[in] fe PETSc basis `FE` object 6817b8c5038SJames Wright @param[out] tabulation PETSc basis tabulation 6827b8c5038SJames Wright @param[out] q_points_1d Quadrature points in libCEED orientation 6837b8c5038SJames Wright @param[out] q_weights_1d Quadrature weights in libCEED orientation 6847b8c5038SJames Wright 6857b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 6867b8c5038SJames Wright **/ 6877b8c5038SJames Wright static inline PetscErrorCode Create1DTabulation_Tensor(MPI_Comm comm, PetscFE fe, PetscTabulation *tabulation, PetscReal **q_points_1d, 6887b8c5038SJames Wright CeedScalar **q_weights_1d) { 6897b8c5038SJames Wright PetscFE fe_1d; 6907b8c5038SJames Wright PetscInt dim, num_comp, Q = -1, q_dim = -1, num_derivatives = 1; 6917b8c5038SJames Wright PetscQuadrature quadrature; 6927b8c5038SJames Wright 6937b8c5038SJames Wright PetscFunctionBeginUser; 6947b8c5038SJames Wright // Create 1D FE 6957b8c5038SJames Wright PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 6967b8c5038SJames Wright PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 6977b8c5038SJames Wright { 6987b8c5038SJames Wright const char *prefix; 6997b8c5038SJames Wright PetscBool is_tensor; 7007b8c5038SJames Wright PetscInt num_comp, order, q_order; 7017b8c5038SJames Wright PetscDualSpace dual_space; 7027b8c5038SJames Wright PetscQuadrature quadrature; 7037b8c5038SJames Wright 7047b8c5038SJames Wright PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix)); 7057b8c5038SJames Wright PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 7067b8c5038SJames Wright PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 7077b8c5038SJames Wright PetscCall(PetscDualSpaceLagrangeGetTensor(dual_space, &is_tensor)); 7087b8c5038SJames Wright PetscCall(PetscDualSpaceGetOrder(dual_space, &order)); 7097b8c5038SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 7107b8c5038SJames Wright PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &q_order, NULL, NULL)); 7117b8c5038SJames Wright PetscCall(PetscFECreateLagrangeFromOptions(comm, 1, num_comp, !is_tensor, order, 7127b8c5038SJames Wright (PetscInt)PetscCeilReal(PetscPowReal(1.0 * q_order, 1.0 / dim)) - 1, prefix, &fe_1d)); 7137b8c5038SJames Wright } 7147b8c5038SJames Wright 7157b8c5038SJames Wright // Get quadrature data 7167b8c5038SJames Wright PetscCall(PetscFEGetQuadrature(fe_1d, &quadrature)); 7177b8c5038SJames Wright PetscCall(GetQuadratureDataP2C(fe_1d, quadrature, &q_dim, NULL, &Q, q_points_1d, q_weights_1d)); 7187b8c5038SJames Wright 7197b8c5038SJames Wright // Compute 1D tabulation 7207b8c5038SJames Wright PetscCall(PetscFECreateTabulation(fe_1d, 1, Q, *q_points_1d, num_derivatives, tabulation)); 7217b8c5038SJames Wright 7227b8c5038SJames Wright // Cleanup 7237b8c5038SJames Wright PetscCall(PetscFEDestroy(&fe_1d)); 7247b8c5038SJames Wright 7257b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7267b8c5038SJames Wright } 7277b8c5038SJames Wright 7287b8c5038SJames Wright /** 7297b8c5038SJames Wright @brief Destroy `CeedBasis` in a `PetscContainer`. 7307b8c5038SJames Wright 7317b8c5038SJames Wright Not collective across MPI processes. 7327b8c5038SJames Wright 7337b8c5038SJames Wright @param[in,out] ctx `CeedBasis` 7347b8c5038SJames Wright 7357b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 7367b8c5038SJames Wright **/ 7377b8c5038SJames Wright static PetscErrorCode DMPlexCeedBasisDestroy(void **ctx) { 7387b8c5038SJames Wright CeedBasis basis = *(CeedBasis *)ctx; 7397b8c5038SJames Wright 7407b8c5038SJames Wright PetscFunctionBegin; 7417b8c5038SJames Wright PetscCallCeed(CeedBasisReturnCeed(basis), CeedBasisDestroy(&basis)); 7427b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7437b8c5038SJames Wright } 7447b8c5038SJames Wright 7457b8c5038SJames Wright /** 7467b8c5038SJames Wright @brief Create `CeedBasis` from `DMPlex`. 7477b8c5038SJames Wright 7487b8c5038SJames Wright The `CeedBasis` is stored in the `DMPlex` object for reuse; 7497b8c5038SJames Wright if the routine is called again with the same arguments, the same `CeedBasis` is returned. 7507b8c5038SJames Wright 7517b8c5038SJames Wright Collective across MPI processes. 7527b8c5038SJames Wright 7537b8c5038SJames Wright @param[in] ceed `Ceed` context 7547b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 7557b8c5038SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 7567b8c5038SJames Wright @param[in] label_value Stratum value 7577b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 7587b8c5038SJames Wright @param[in] dm_field Index of `DMPlex` field 7597b8c5038SJames Wright @param[out] basis `CeedBasis` for `DMPlex` 7607b8c5038SJames Wright 7617b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 7627b8c5038SJames Wright **/ 7637b8c5038SJames Wright PetscErrorCode DMPlexCeedBasisCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field, 7647b8c5038SJames Wright CeedBasis *basis) { 7657b8c5038SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 7667b8c5038SJames Wright char container_name[1024]; 7677b8c5038SJames Wright CeedBasis container_basis; 7687b8c5038SJames Wright 7697b8c5038SJames Wright PetscFunctionBeginUser; 7707b8c5038SJames Wright { 7717b8c5038SJames Wright const char *label_name = NULL; 7727b8c5038SJames Wright 7737b8c5038SJames Wright if (domain_label) PetscCall(PetscObjectGetName((PetscObject)domain_label, &label_name)); 7747b8c5038SJames Wright PetscCall(PetscSNPrintf(container_name, sizeof(container_name), 7757b8c5038SJames Wright "DM CeedBasis - label: %s label value: %" PetscInt_FMT " height: %" PetscInt_FMT " DM field: %" PetscInt_FMT, 7767b8c5038SJames Wright label_name ? label_name : "NONE", label_value, height, dm_field)); 7777b8c5038SJames Wright } 7787b8c5038SJames Wright PetscCall(PetscObjectContainerQuery((PetscObject)dm, container_name, (void **)&container_basis)); 7797b8c5038SJames Wright 7807b8c5038SJames Wright if (container_basis) { 7817b8c5038SJames Wright // Copy existing basis 7827b8c5038SJames Wright *basis = NULL; 7837b8c5038SJames Wright PetscCallCeed(ceed, CeedBasisReferenceCopy(container_basis, basis)); 7847b8c5038SJames Wright } else { 7857b8c5038SJames Wright PetscDS ds; 7867b8c5038SJames Wright PetscFE fe; 7877b8c5038SJames Wright PetscDualSpace dual_space; 7887b8c5038SJames Wright PetscQuadrature quadrature; 7897b8c5038SJames Wright PetscBool is_tensor = PETSC_TRUE; 7907b8c5038SJames Wright PetscInt ds_field = -1; 7917b8c5038SJames Wright 7927b8c5038SJames Wright // Get element information 7937b8c5038SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 7947b8c5038SJames Wright PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 7957b8c5038SJames Wright PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 7967b8c5038SJames Wright PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 7977b8c5038SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 7987b8c5038SJames Wright PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 7997b8c5038SJames Wright PetscCall(PetscDualSpaceLagrangeGetTensor(dual_space, &is_tensor)); 8007b8c5038SJames Wright 8017b8c5038SJames Wright // Build libCEED basis 8027b8c5038SJames Wright PetscBool use_nontensor = PETSC_FALSE; 8037b8c5038SJames Wright 8047b8c5038SJames Wright PetscOptionsBegin(comm, NULL, "Tensor Basis as Nontensor Check", NULL); 8057b8c5038SJames Wright PetscCall(PetscOptionsBool("-ceed_basis_as_nontensor", "", NULL, use_nontensor, &use_nontensor, NULL)); 8067b8c5038SJames Wright PetscOptionsEnd(); 8077b8c5038SJames Wright 8087b8c5038SJames Wright if (!is_tensor || use_nontensor) { 8097b8c5038SJames Wright PetscTabulation basis_tabulation; 8107b8c5038SJames Wright PetscInt num_derivatives = 1, face = 0; 8117b8c5038SJames Wright CeedScalar *q_points, *q_weights, *interp, *grad; 8127b8c5038SJames Wright CeedElemTopology elem_topo; 8137b8c5038SJames Wright 8147b8c5038SJames Wright { 8157b8c5038SJames Wright DMPolytopeType cell_type; 8167b8c5038SJames Wright 8177b8c5038SJames Wright PetscCall(GetGlobalDMPlexPolytopeType(dm, domain_label, label_value, height, &cell_type)); 8187b8c5038SJames Wright elem_topo = PolytopeTypePetscToCeed(cell_type); 8197b8c5038SJames Wright PetscCheck(elem_topo, comm, PETSC_ERR_SUP, "DMPlex topology not supported: %s", DMPolytopeTypes[cell_type]); 8207b8c5038SJames Wright } 8217b8c5038SJames Wright 8227b8c5038SJames Wright // Compute basis tabulation 8237b8c5038SJames Wright PetscCall(GetQuadratureDataP2C(fe, quadrature, NULL, NULL, NULL, &q_points, &q_weights)); 8247b8c5038SJames Wright PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation)); 8257b8c5038SJames Wright PetscCall(ComputeFieldTabulationP2C(dm, dm_field, face, basis_tabulation, &interp, &grad)); 8267b8c5038SJames Wright { 8277b8c5038SJames Wright PetscInt num_comp = basis_tabulation->Nc, P = basis_tabulation->Nb / num_comp, Q = basis_tabulation->Np; 8287b8c5038SJames Wright 8297b8c5038SJames Wright PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis)); 8307b8c5038SJames Wright } 8317b8c5038SJames Wright 8327b8c5038SJames Wright PetscCall(PetscFree(q_points)); 8337b8c5038SJames Wright PetscCall(PetscFree(q_weights)); 8347b8c5038SJames Wright PetscCall(PetscFree(interp)); 8357b8c5038SJames Wright PetscCall(PetscFree(grad)); 8367b8c5038SJames Wright } else { 8377b8c5038SJames Wright PetscInt P_1d, Q_1d, num_comp, dim; 8387b8c5038SJames Wright PetscTabulation basis_tabulation; 8397b8c5038SJames Wright CeedScalar *q_points_1d, *q_weights_1d, *interp_1d, *grad_1d; 8407b8c5038SJames Wright 8417b8c5038SJames Wright PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 8427b8c5038SJames Wright PetscCall(Create1DTabulation_Tensor(comm, fe, &basis_tabulation, &q_points_1d, &q_weights_1d)); 8437b8c5038SJames Wright num_comp = basis_tabulation->Nc; 8447b8c5038SJames Wright P_1d = basis_tabulation->Nb / num_comp; 8457b8c5038SJames Wright Q_1d = basis_tabulation->Np; 8467b8c5038SJames Wright PetscCall(ComputeFieldTabulationP2C(dm, dm_field, 0, basis_tabulation, &interp_1d, &grad_1d)); 8477b8c5038SJames Wright PetscCallCeed(ceed, CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_points_1d, q_weights_1d, basis)); 8487b8c5038SJames Wright 8497b8c5038SJames Wright // Cleanup 8507b8c5038SJames Wright PetscCall(PetscFree(q_points_1d)); 8517b8c5038SJames Wright PetscCall(PetscFree(q_weights_1d)); 8527b8c5038SJames Wright PetscCall(PetscFree(interp_1d)); 8537b8c5038SJames Wright PetscCall(PetscFree(grad_1d)); 8547b8c5038SJames Wright PetscCall(PetscTabulationDestroy(&basis_tabulation)); 8557b8c5038SJames Wright } 8567b8c5038SJames Wright // Set in container 8577b8c5038SJames Wright PetscCallCeed(ceed, CeedBasisReferenceCopy(*basis, &container_basis)); 8587b8c5038SJames Wright PetscCall(PetscObjectContainerCompose((PetscObject)dm, container_name, container_basis, DMPlexCeedBasisDestroy)); 8597b8c5038SJames Wright } 8607b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 8617b8c5038SJames Wright } 8627b8c5038SJames Wright 8637b8c5038SJames Wright /** 8647b8c5038SJames Wright @brief Create `CeedBasis` for coordinate space of `DMPlex`. 8657b8c5038SJames Wright 8667b8c5038SJames Wright The `CeedBasis` is stored in the coordinate `DMPlex` object for reuse; 8677b8c5038SJames Wright if the routine is called again with the same arguments, the same `CeedBasis` is returned. 8687b8c5038SJames Wright 8697b8c5038SJames Wright Collective across MPI processes. 8707b8c5038SJames Wright 8717b8c5038SJames Wright @param[in] ceed `Ceed` context 8727b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 8737b8c5038SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 8747b8c5038SJames Wright @param[in] label_value Stratum value 8757b8c5038SJames Wright @param[in] height Height of `DMPlex` topology 8767b8c5038SJames Wright @param[out] basis `CeedBasis` for coordinate space of `DMPlex` 8777b8c5038SJames Wright **/ 8787b8c5038SJames Wright PetscErrorCode DMPlexCeedBasisCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, CeedBasis *basis) { 8797b8c5038SJames Wright DM dm_coord = NULL; 8807b8c5038SJames Wright 8817b8c5038SJames Wright PetscFunctionBeginUser; 882*77efe802SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 8837b8c5038SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, dm_coord, domain_label, label_value, height, 0, basis)); 8847b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 8857b8c5038SJames Wright } 8867b8c5038SJames Wright 887*77efe802SJames Wright typedef struct CeedVectorAndState_ *CeedVectorAndState; 888*77efe802SJames Wright struct CeedVectorAndState_ { 889*77efe802SJames Wright CeedVector vector; // !< `CeedVector` containing copy of data from a `Vec` 890*77efe802SJames Wright PetscObjectState state; // !< Object state of the `Vec` when the data was copied 891*77efe802SJames Wright }; 892*77efe802SJames Wright 893*77efe802SJames Wright /** 894*77efe802SJames Wright @brief Destroy `CeedVectorAndState` in a `PetscContainer`. 895*77efe802SJames Wright 896*77efe802SJames Wright Not collective across MPI processes. 897*77efe802SJames Wright 898*77efe802SJames Wright @param[in,out] ctx `CeedVectorAndState` 899*77efe802SJames Wright 900*77efe802SJames Wright @return An error code: 0 - success, otherwise - failure 901*77efe802SJames Wright **/ 902*77efe802SJames Wright static PetscErrorCode CeedVectorAndStateDestroy(void **ctx) { 903*77efe802SJames Wright CeedVectorAndState vec_and_state = *(CeedVectorAndState *)ctx; 904*77efe802SJames Wright 905*77efe802SJames Wright PetscFunctionBegin; 906*77efe802SJames Wright PetscCallCeed(CeedVectorReturnCeed(vec_and_state->vector), CeedVectorDestroy(&vec_and_state->vector)); 907*77efe802SJames Wright PetscCall(PetscFree(vec_and_state)); 908*77efe802SJames Wright *ctx = NULL; 909*77efe802SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 910*77efe802SJames Wright } 911*77efe802SJames Wright 912*77efe802SJames Wright /** 913*77efe802SJames Wright @brief Get Ceed objects for coordinate field of DM 914*77efe802SJames Wright 915*77efe802SJames Wright @param[in] ceed `Ceed` context 916*77efe802SJames Wright @param[in] dm `DMPlex` holding mesh 917*77efe802SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 918*77efe802SJames Wright @param[in] label_value Stratum value 919*77efe802SJames Wright @param[in] height Height of `DMPlex` topology 920*77efe802SJames Wright @param[out] elem_restr `CeedElemRestriction` for coodinates of `DMPlex`, or `NULL` 921*77efe802SJames Wright @param[out] basis `CeedBasis` for coordinate space of `DMPlex`, or `NULL` 922*77efe802SJames Wright @param[out] vector `CeedVector` for coordinates of `DMPlex`, or `NULL` 923*77efe802SJames Wright **/ 924*77efe802SJames Wright PetscErrorCode DMPlexCeedCoordinateCreateField(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, 925*77efe802SJames Wright CeedElemRestriction *elem_restr, CeedBasis *basis, CeedVector *vector) { 926*77efe802SJames Wright CeedElemRestriction elem_restr_ = NULL; 927*77efe802SJames Wright 928*77efe802SJames Wright PetscFunctionBeginUser; 929*77efe802SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_)); 930*77efe802SJames Wright if (elem_restr) { 931*77efe802SJames Wright *elem_restr = NULL; 932*77efe802SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(elem_restr_, elem_restr)); 933*77efe802SJames Wright } 934*77efe802SJames Wright if (basis) PetscCall(DMPlexCeedBasisCoordinateCreate(ceed, dm, domain_label, label_value, height, basis)); 935*77efe802SJames Wright if (vector) { 936*77efe802SJames Wright DM cdm; 937*77efe802SJames Wright char container_name[] = "DM Coordinate CeedVectorAndState"; 938*77efe802SJames Wright CeedVectorAndState vec_and_state; 939*77efe802SJames Wright Vec X_loc; 940*77efe802SJames Wright PetscObjectState X_loc_state; 941*77efe802SJames Wright 942*77efe802SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 943*77efe802SJames Wright if (cdm) { 944*77efe802SJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 945*77efe802SJames Wright } else { 946*77efe802SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 947*77efe802SJames Wright } 948*77efe802SJames Wright PetscCall(VecGetState(X_loc, &X_loc_state)); 949*77efe802SJames Wright PetscCall(PetscObjectContainerQuery((PetscObject)X_loc, container_name, (void **)&vec_and_state)); 950*77efe802SJames Wright 951*77efe802SJames Wright if (vec_and_state && vec_and_state->state == X_loc_state) { 952*77efe802SJames Wright // Copy existing vector 953*77efe802SJames Wright *vector = NULL; 954*77efe802SJames Wright PetscCallCeed(ceed, CeedVectorReferenceCopy(vec_and_state->vector, vector)); 955*77efe802SJames Wright } else { 956*77efe802SJames Wright PetscCall(CeedElemRestrictionCreateVector(elem_restr_, vector, NULL)); 957*77efe802SJames Wright PetscCall(VecCopyPetscToCeed(X_loc, *vector)); 958*77efe802SJames Wright 959*77efe802SJames Wright // Set in container 960*77efe802SJames Wright PetscCall(PetscNew(&vec_and_state)); 961*77efe802SJames Wright PetscCallCeed(ceed, CeedVectorReferenceCopy(*vector, &vec_and_state->vector)); 962*77efe802SJames Wright vec_and_state->state = X_loc_state; 963*77efe802SJames Wright PetscCall(PetscObjectContainerCompose((PetscObject)X_loc, container_name, vec_and_state, CeedVectorAndStateDestroy)); 964*77efe802SJames Wright } 965*77efe802SJames Wright } 966*77efe802SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_)); 967*77efe802SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 968*77efe802SJames Wright } 969*77efe802SJames Wright 9707b8c5038SJames Wright /** 9717b8c5038SJames Wright @brief Create `CeedBasis` for cell to face quadrature space evaluation from `DMPlex`. 9727b8c5038SJames Wright 9737b8c5038SJames Wright Collective across MPI processes. 9747b8c5038SJames Wright 9757b8c5038SJames Wright @param[in] ceed `Ceed` context 9767b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 9777b8c5038SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 9787b8c5038SJames Wright @param[in] label_value Stratum value 9797b8c5038SJames Wright @param[in] face Index of face 9807b8c5038SJames Wright @param[in] dm_field Index of `DMPlex` field 9817b8c5038SJames Wright @param[out] basis `CeedBasis` for cell to face evaluation 9827b8c5038SJames Wright **/ 9837b8c5038SJames Wright PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field, 9847b8c5038SJames Wright CeedBasis *basis) { 9857b8c5038SJames Wright PetscInt ds_field = -1, height = 0; 9867b8c5038SJames Wright PetscDS ds; 9877b8c5038SJames Wright PetscFE fe; 9887b8c5038SJames Wright PetscQuadrature face_quadrature; 9897b8c5038SJames Wright 9907b8c5038SJames Wright PetscFunctionBeginUser; 9917b8c5038SJames Wright // Get element information 9927b8c5038SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 9937b8c5038SJames Wright PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 9947b8c5038SJames Wright PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 9957b8c5038SJames Wright PetscCall(PetscFEGetFaceQuadrature(fe, &face_quadrature)); 9967b8c5038SJames Wright 9977b8c5038SJames Wright { // Build libCEED basis 9987b8c5038SJames Wright PetscInt num_derivatives = 1, num_comp, P, Q = -1; 9997b8c5038SJames Wright PetscTabulation basis_tabulation; 10007b8c5038SJames Wright CeedScalar *q_points, *q_weights, *interp, *grad; 10017b8c5038SJames Wright CeedElemTopology elem_topo; 10027b8c5038SJames Wright 10037b8c5038SJames Wright { // Verify compatible element topology 10047b8c5038SJames Wright DMPolytopeType cell_type; 10057b8c5038SJames Wright 10067b8c5038SJames Wright PetscCall(GetGlobalDMPlexPolytopeType(dm, domain_label, label_value, height, &cell_type)); 10077b8c5038SJames Wright elem_topo = PolytopeTypePetscToCeed(cell_type); 10087b8c5038SJames Wright PetscCheck(elem_topo, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported: %s", DMPolytopeTypes[cell_type]); 10097b8c5038SJames Wright } 10107b8c5038SJames Wright 10117b8c5038SJames Wright // Compute basis tabulation 10127b8c5038SJames Wright PetscCall(GetQuadratureDataP2C(fe, face_quadrature, NULL, NULL, &Q, &q_points, &q_weights)); 10137b8c5038SJames Wright PetscCall(PetscFEGetFaceTabulation(fe, num_derivatives, &basis_tabulation)); 10147b8c5038SJames Wright num_comp = basis_tabulation->Nc; 10157b8c5038SJames Wright P = basis_tabulation->Nb / num_comp; 10167b8c5038SJames Wright PetscCall(ComputeFieldTabulationP2C(dm, dm_field, face, basis_tabulation, &interp, &grad)); 10177b8c5038SJames Wright PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis)); 10187b8c5038SJames Wright 10197b8c5038SJames Wright PetscCall(PetscFree(q_points)); 10207b8c5038SJames Wright PetscCall(PetscFree(q_weights)); 10217b8c5038SJames Wright PetscCall(PetscFree(interp)); 10227b8c5038SJames Wright PetscCall(PetscFree(grad)); 10237b8c5038SJames Wright } 10247b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 10257b8c5038SJames Wright } 10267b8c5038SJames Wright 10277b8c5038SJames Wright /** 10287b8c5038SJames Wright @brief Create `CeedBasis` for cell to face quadrature space evaluation on coordinate space of `DMPlex`. 10297b8c5038SJames Wright 10307b8c5038SJames Wright Collective across MPI processes. 10317b8c5038SJames Wright 10327b8c5038SJames Wright @param[in] ceed `Ceed` context 10337b8c5038SJames Wright @param[in] dm `DMPlex` holding mesh 10347b8c5038SJames Wright @param[in] domain_label `DMLabel` for `DMPlex` domain 10357b8c5038SJames Wright @param[in] label_value Stratum value 10367b8c5038SJames Wright @param[in] face Index of face 10377b8c5038SJames Wright @param[out] basis `CeedBasis` for cell to face evaluation on coordinate space of `DMPlex` 10387b8c5038SJames Wright **/ 10397b8c5038SJames Wright PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, 10407b8c5038SJames Wright CeedBasis *basis) { 10417b8c5038SJames Wright DM dm_coord; 10427b8c5038SJames Wright 10437b8c5038SJames Wright PetscFunctionBeginUser; 10447b8c5038SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 10457b8c5038SJames Wright if (!dm_coord) { 10467b8c5038SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 10477b8c5038SJames Wright } 10487b8c5038SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, dm_coord, domain_label, label_value, face, 0, basis)); 10497b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 10507b8c5038SJames Wright } 10517b8c5038SJames Wright 10527b8c5038SJames Wright /** 10537b8c5038SJames Wright @brief Setup `DM` with FE space of appropriate degree 10547b8c5038SJames Wright 10557b8c5038SJames Wright Must be followed by `DMSetupByOrderEnd_FEM` 10567b8c5038SJames Wright 10577b8c5038SJames Wright @param[in] setup_faces Flag to setup face geometry 10587b8c5038SJames Wright @param[in] setup_coords Flag to setup coordinate spaces 10597b8c5038SJames Wright @param[in] degree Polynomial orders of field 10607b8c5038SJames Wright @param[in] coord_order Polynomial order of coordinate basis, or `PETSC_DECIDE` for default 10617b8c5038SJames Wright @param[in] q_extra Additional quadrature order 10627b8c5038SJames Wright @param[in] num_fields Number of fields in solution vector 10637b8c5038SJames Wright @param[in] field_sizes Array of number of components for each field 10647b8c5038SJames Wright @param[out] dm `DM` to setup 10657b8c5038SJames Wright 10667b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 10677b8c5038SJames Wright **/ 10687b8c5038SJames Wright PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 10697b8c5038SJames Wright PetscInt num_fields, const PetscInt *field_sizes, DM dm) { 10707b8c5038SJames Wright PetscInt dim, q_order = degree + q_extra; 10717b8c5038SJames Wright PetscBool is_simplex = PETSC_TRUE; 10727b8c5038SJames Wright PetscFE fe; 10737b8c5038SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 10747b8c5038SJames Wright 10757b8c5038SJames Wright PetscFunctionBeginUser; 10767b8c5038SJames Wright PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 10777b8c5038SJames Wright 10787b8c5038SJames Wright // Setup DM 10797b8c5038SJames Wright PetscCall(DMGetDimension(dm, &dim)); 10807b8c5038SJames Wright for (PetscInt i = 0; i < num_fields; i++) { 10817b8c5038SJames Wright PetscFE fe_face; 10827b8c5038SJames Wright PetscInt q_order = degree + q_extra; 10837b8c5038SJames Wright 10847b8c5038SJames Wright PetscCall(PetscFECreateLagrange(comm, dim, field_sizes[i], is_simplex, degree, q_order, &fe)); 10857b8c5038SJames Wright if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe_face)); 10867b8c5038SJames Wright PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 10877b8c5038SJames Wright PetscCall(PetscFEDestroy(&fe)); 10887b8c5038SJames Wright } 10897b8c5038SJames Wright PetscCall(DMCreateDS(dm)); 10907b8c5038SJames Wright 10917b8c5038SJames Wright // Project coordinates to enrich quadrature space 10927b8c5038SJames Wright if (setup_coords) { 10937b8c5038SJames Wright DM dm_coord; 10947b8c5038SJames Wright PetscDS ds_coord; 10957b8c5038SJames Wright PetscFE fe_coord_current, fe_coord_new, fe_coord_face_new; 10967b8c5038SJames Wright PetscDualSpace fe_coord_dual_space; 10977b8c5038SJames Wright PetscInt fe_coord_order, num_comp_coord; 10987b8c5038SJames Wright 10997b8c5038SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 11007b8c5038SJames Wright PetscCall(DMGetCoordinateDim(dm, &num_comp_coord)); 11017b8c5038SJames Wright PetscCall(DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord, NULL)); 11027b8c5038SJames Wright PetscCall(PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current)); 11037b8c5038SJames Wright PetscCall(PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space)); 11047b8c5038SJames Wright PetscCall(PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order)); 11057b8c5038SJames Wright 11067b8c5038SJames Wright // Create FE for coordinates 11077b8c5038SJames Wright if (coord_order != PETSC_DECIDE) fe_coord_order = coord_order; 11087b8c5038SJames Wright PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, q_order, &fe_coord_new)); 11097b8c5038SJames Wright if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe_coord_new, 1, &fe_coord_face_new)); 11107b8c5038SJames Wright PetscCall(DMSetCoordinateDisc(dm, fe_coord_new, PETSC_FALSE, PETSC_TRUE)); 11117b8c5038SJames Wright PetscCall(PetscFEDestroy(&fe_coord_new)); 11127b8c5038SJames Wright } 11137b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 11147b8c5038SJames Wright } 11157b8c5038SJames Wright 11167b8c5038SJames Wright /** 11177b8c5038SJames Wright @brief Finish setting up `DM` 11187b8c5038SJames Wright 11197b8c5038SJames Wright Must be called after `DMSetupByOrderBegin_FEM` and all strong BCs have be declared via `DMAddBoundaries`. 11207b8c5038SJames Wright 11217b8c5038SJames Wright @param[in] setup_coords Flag to setup coordinate spaces 11227b8c5038SJames Wright @param[out] dm `DM` to setup 11237b8c5038SJames Wright 11247b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 11257b8c5038SJames Wright **/ 11267b8c5038SJames Wright PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm) { 11277b8c5038SJames Wright PetscBool is_simplex; 11287b8c5038SJames Wright 11297b8c5038SJames Wright PetscFunctionBeginUser; 11307b8c5038SJames Wright PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 11317b8c5038SJames Wright // Set tensor permutation if needed 11327b8c5038SJames Wright if (!is_simplex) { 11337b8c5038SJames Wright PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 11347b8c5038SJames Wright if (setup_coords) { 11357b8c5038SJames Wright DM dm_coord; 11367b8c5038SJames Wright 11377b8c5038SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 11387b8c5038SJames Wright PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 11397b8c5038SJames Wright } 11407b8c5038SJames Wright } 11417b8c5038SJames Wright PetscCall(DMLocalizeCoordinates(dm)); // Must localize after tensor closure setting 11427b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 11437b8c5038SJames Wright } 11447b8c5038SJames Wright 11457b8c5038SJames Wright /** 11467b8c5038SJames Wright @brief Setup `DM` with FE space of appropriate degree with no boundary conditions 11477b8c5038SJames Wright 11487b8c5038SJames Wright Calls `DMSetupByOrderBegin_FEM` and `DMSetupByOrderEnd_FEM` successively 11497b8c5038SJames Wright 11507b8c5038SJames Wright @param[in] setup_faces Flag to setup face geometry 11517b8c5038SJames Wright @param[in] setup_coords Flag to setup coordinate spaces 11527b8c5038SJames Wright @param[in] degree Polynomial orders of field 11537b8c5038SJames Wright @param[in] coord_order Polynomial order of coordinate basis, or `PETSC_DECIDE` for default 11547b8c5038SJames Wright @param[in] q_extra Additional quadrature order 11557b8c5038SJames Wright @param[in] num_fields Number of fields in solution vector 11567b8c5038SJames Wright @param[in] field_sizes Array of number of components for each field 11577b8c5038SJames Wright @param[out] dm `DM` to setup 11587b8c5038SJames Wright 11597b8c5038SJames Wright @return An error code: 0 - success, otherwise - failure 11607b8c5038SJames Wright **/ 11617b8c5038SJames Wright PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra, 11627b8c5038SJames Wright PetscInt num_fields, const PetscInt *field_sizes, DM dm) { 11637b8c5038SJames Wright PetscFunctionBeginUser; 11647b8c5038SJames Wright PetscCall(DMSetupByOrderBegin_FEM(setup_faces, setup_coords, degree, coord_order, q_extra, num_fields, field_sizes, dm)); 11657b8c5038SJames Wright PetscCall(DMSetupByOrderEnd_FEM(setup_coords, dm)); 11667b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 11677b8c5038SJames Wright } 11687b8c5038SJames Wright 11697b8c5038SJames Wright /** 11707b8c5038SJames Wright @brief Get the points in a label stratum which have requested height 11717b8c5038SJames Wright 11727b8c5038SJames Wright 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: 11737b8c5038SJames Wright ``` 11747b8c5038SJames Wright PetscCall(DMGetStratumISAtHeight(dm, "Face Sets", 3, 1, &is_face_points)); 11757b8c5038SJames Wright ``` 11767b8c5038SJames Wright 11777b8c5038SJames Wright @param[in] dm DM which has the label 11787b8c5038SJames Wright @param[in] name Name fo the label 11797b8c5038SJames Wright @param[in] value Label value of points to return 11807b8c5038SJames Wright @param[in] height Height of points to return 11817b8c5038SJames Wright @param[out] points `IS` of points, or `NULL` is there are no points in that stratum at height 11827b8c5038SJames Wright **/ 11837b8c5038SJames Wright static PetscErrorCode DMGetStratumISAtHeight(DM dm, const char name[], PetscInt value, PetscInt height, IS *points) { 11847b8c5038SJames Wright PetscInt depth; 11857b8c5038SJames Wright DMLabel depth_label; 11867b8c5038SJames Wright IS label_is, depth_is; 11877b8c5038SJames Wright 11887b8c5038SJames Wright PetscFunctionBeginUser; 11897b8c5038SJames Wright PetscCall(DMPlexGetDepth(dm, &depth)); 11907b8c5038SJames Wright PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 11917b8c5038SJames Wright PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is)); 11927b8c5038SJames Wright PetscCall(DMGetStratumIS(dm, name, value, &label_is)); 11937b8c5038SJames Wright if (label_is == NULL || depth_is == NULL) *points = NULL; 11947b8c5038SJames Wright else PetscCall(ISIntersect(depth_is, label_is, points)); 11957b8c5038SJames Wright PetscCall(ISDestroy(&depth_is)); 11967b8c5038SJames Wright PetscCall(ISDestroy(&label_is)); 11977b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 11987b8c5038SJames Wright } 11997b8c5038SJames Wright 12007b8c5038SJames Wright /** 12017b8c5038SJames Wright @brief Create a label with separate values for the `FE` orientations for all cells with this material for the `DMPlex` face. 12027b8c5038SJames Wright 12037b8c5038SJames Wright Collective across MPI processes. 12047b8c5038SJames Wright 12057b8c5038SJames Wright @param[in,out] dm `DMPlex` holding mesh 12067b8c5038SJames Wright @param[in] dm_face Face number on `DMPlex` 12077b8c5038SJames Wright @param[out] face_label_name Label name for newly created label, or `NULL` if no elements match the `material` and `dm_face` 12087b8c5038SJames Wright **/ 12097b8c5038SJames Wright PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, char **face_label_name) { 12107b8c5038SJames Wright PetscInt num_face_points = 0, face_height = 1; 12117b8c5038SJames Wright const PetscInt *face_points; 12127b8c5038SJames Wright IS is_face_points; 12137b8c5038SJames Wright DMLabel face_label = NULL; 12147b8c5038SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 12157b8c5038SJames Wright 12167b8c5038SJames Wright PetscFunctionBeginUser; 12177b8c5038SJames Wright { // Create label name and label 12187b8c5038SJames Wright PetscSizeT label_name_len = 64; 12197b8c5038SJames Wright PetscBool has_label_already; 12207b8c5038SJames Wright 12217b8c5038SJames Wright // Face label 12227b8c5038SJames Wright PetscCall(PetscCalloc1(label_name_len, face_label_name)); 12237b8c5038SJames Wright PetscCall(PetscSNPrintf(*face_label_name, label_name_len, "Orientation for DM face %" PetscInt_FMT, dm_face)); 12247b8c5038SJames Wright PetscCall(DMHasLabel(dm, *face_label_name, &has_label_already)); 12257b8c5038SJames Wright if (has_label_already) PetscFunctionReturn(PETSC_SUCCESS); 12267b8c5038SJames Wright 12277b8c5038SJames Wright PetscCall(DMCreateLabel(dm, *face_label_name)); 12287b8c5038SJames Wright PetscCall(DMGetLabel(dm, *face_label_name, &face_label)); 12297b8c5038SJames Wright } 12307b8c5038SJames Wright 12317b8c5038SJames Wright PetscCall(DMGetStratumISAtHeight(dm, "Face Sets", dm_face, face_height, &is_face_points)); 12327b8c5038SJames Wright if (is_face_points) { 12337b8c5038SJames Wright PetscCall(ISGetLocalSize(is_face_points, &num_face_points)); 12347b8c5038SJames Wright PetscCall(ISGetIndices(is_face_points, &face_points)); 12357b8c5038SJames Wright } 12367b8c5038SJames Wright 12377b8c5038SJames Wright for (PetscInt face_index = 0; face_index < num_face_points; face_index++) { 12387b8c5038SJames Wright const PetscInt *face_support, *cell_cone; 12397b8c5038SJames Wright PetscInt face_support_size = 0, cell_cone_size = 0, fe_face_on_cell = -1; 12407b8c5038SJames Wright PetscInt face_point = face_points[face_index]; 12417b8c5038SJames Wright 12427b8c5038SJames Wright // Get cell supporting face 12437b8c5038SJames Wright PetscCall(DMPlexGetSupport(dm, face_point, &face_support)); 12447b8c5038SJames Wright PetscCall(DMPlexGetSupportSize(dm, face_point, &face_support_size)); 12457b8c5038SJames Wright PetscCheck(face_support_size == 1, comm, PETSC_ERR_SUP, "Expected one cell in support of exterior face, but got %" PetscInt_FMT " cells", 12467b8c5038SJames Wright face_support_size); 12477b8c5038SJames Wright PetscInt cell_point = face_support[0]; 12487b8c5038SJames Wright 12497b8c5038SJames Wright // Find face number 12507b8c5038SJames Wright PetscCall(DMPlexGetCone(dm, cell_point, &cell_cone)); 12517b8c5038SJames Wright PetscCall(DMPlexGetConeSize(dm, cell_point, &cell_cone_size)); 12527b8c5038SJames Wright for (PetscInt i = 0; i < cell_cone_size; i++) { 12537b8c5038SJames Wright if (cell_cone[i] == face_point) fe_face_on_cell = i; 12547b8c5038SJames Wright } 12557b8c5038SJames Wright PetscCheck(fe_face_on_cell != -1, comm, PETSC_ERR_SUP, "Could not find face %" PetscInt_FMT " in cone of its support", face_point); 12567b8c5038SJames Wright 12577b8c5038SJames Wright // Add to label 12587b8c5038SJames Wright PetscCall(DMLabelSetValue(face_label, face_point, fe_face_on_cell)); 12597b8c5038SJames Wright } 12607b8c5038SJames Wright PetscCall(DMPlexLabelAddFaceCells(dm, face_label)); 12617b8c5038SJames Wright PetscCall(ISDestroy(&is_face_points)); 12627b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 12637b8c5038SJames Wright } 12647b8c5038SJames Wright 12657b8c5038SJames Wright /** 12667b8c5038SJames Wright @brief Creates an array of all the values set for a `DMLabel` 12677b8c5038SJames Wright 12687b8c5038SJames Wright Because `DMLabelGetValueIS` returns label values locally, not globally. 12697b8c5038SJames Wright 12707b8c5038SJames Wright Collective across MPI processes. 12717b8c5038SJames Wright 12727b8c5038SJames Wright @param[in] dm `DM` with the label 12737b8c5038SJames Wright @param[in] label `DMLabel` to get values of 12747b8c5038SJames Wright @param[out] num_values Total number of values 12757b8c5038SJames Wright @param[out] value_array Array of label values, must be freed by user 12767b8c5038SJames Wright **/ 12777b8c5038SJames Wright PetscErrorCode DMLabelCreateGlobalValueArray(DM dm, DMLabel label, PetscInt *num_values, PetscInt **value_array) { 12787b8c5038SJames Wright PetscInt num_values_local, minmax_values[2], minmax_values_loc[2]; 12797b8c5038SJames Wright PetscInt *values_local = NULL; 12807b8c5038SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 12817b8c5038SJames Wright PetscSegBuffer seg; 12827b8c5038SJames Wright PetscCount seg_size; 12837b8c5038SJames Wright 12847b8c5038SJames Wright PetscFunctionBeginUser; 12857b8c5038SJames Wright { // Extract array of local DMLabel values so it can be sorted 12867b8c5038SJames Wright IS is_values; 12877b8c5038SJames Wright const PetscInt *is_values_local = NULL; 12887b8c5038SJames Wright 12897b8c5038SJames Wright PetscCall(DMLabelGetValueIS(label, &is_values)); 12907b8c5038SJames Wright PetscCall(ISGetLocalSize(is_values, &num_values_local)); 12917b8c5038SJames Wright 12927b8c5038SJames Wright PetscCall(ISGetIndices(is_values, &is_values_local)); 12937b8c5038SJames Wright PetscCall(PetscMalloc1(num_values_local, &values_local)); 12947b8c5038SJames Wright PetscCall(PetscArraycpy(values_local, is_values_local, num_values_local)); 12957b8c5038SJames Wright PetscCall(ISRestoreIndices(is_values, &is_values_local)); 12967b8c5038SJames Wright PetscCall(ISDestroy(&is_values)); 12977b8c5038SJames Wright } 12987b8c5038SJames Wright 12997b8c5038SJames Wright if (num_values_local) { 13007b8c5038SJames Wright PetscCall(PetscSortInt(num_values_local, values_local)); 13017b8c5038SJames Wright minmax_values_loc[0] = values_local[0]; 13027b8c5038SJames Wright minmax_values_loc[1] = values_local[num_values_local - 1]; 13037b8c5038SJames Wright } else { 13047b8c5038SJames Wright // Rank has no set label values, therefore set local min/max to have no effect on global min/max 13057b8c5038SJames Wright minmax_values_loc[0] = PETSC_INT_MAX; 13067b8c5038SJames Wright minmax_values_loc[1] = PETSC_INT_MIN; 13077b8c5038SJames Wright } 13087b8c5038SJames Wright PetscCall(PetscGlobalMinMaxInt(comm, minmax_values_loc, minmax_values)); 13097b8c5038SJames Wright 13107b8c5038SJames Wright PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 16, &seg)); 13117b8c5038SJames Wright for (PetscInt value = minmax_values[0]; value <= minmax_values[1]; value++) { 13127b8c5038SJames Wright PetscInt location_local = -1, location_global = -1; 13137b8c5038SJames Wright 13147b8c5038SJames Wright PetscCall(PetscFindInt(value, num_values_local, values_local, &location_local)); 13157b8c5038SJames Wright PetscCallMPI(MPIU_Allreduce(&location_local, &location_global, 1, MPIU_INT, MPI_MAX, comm)); 13167b8c5038SJames Wright if (location_global < 0) continue; 13177b8c5038SJames Wright else { 13187b8c5038SJames Wright PetscInt *seg_slot; 13197b8c5038SJames Wright 13207b8c5038SJames Wright PetscCall(PetscSegBufferGetInts(seg, 1, &seg_slot)); 13217b8c5038SJames Wright *seg_slot = value; 13227b8c5038SJames Wright } 13237b8c5038SJames Wright } 13247b8c5038SJames Wright PetscCall(PetscSegBufferGetSize(seg, &seg_size)); 13257b8c5038SJames Wright *num_values = (PetscInt)seg_size; 13267b8c5038SJames Wright PetscCall(PetscSegBufferExtractAlloc(seg, value_array)); 13277b8c5038SJames Wright PetscCall(PetscSegBufferDestroy(&seg)); 13287b8c5038SJames Wright if (values_local) PetscCall(PetscFree(values_local)); 13297b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 13307b8c5038SJames Wright } 13317b8c5038SJames Wright 13327b8c5038SJames Wright /** 13337b8c5038SJames Wright @brief Get the number of components for `dm`'s coordinate field 13347b8c5038SJames Wright 13357b8c5038SJames Wright @param[in] dm DM 13367b8c5038SJames Wright @param[out] num_comp Number of components for the coordinate field 13377b8c5038SJames Wright **/ 13387b8c5038SJames Wright PetscErrorCode DMGetCoordinateNumComps(DM dm, PetscInt *num_comp) { 13397b8c5038SJames Wright DM dm_coord; 13407b8c5038SJames Wright PetscSection section_coord; 13417b8c5038SJames Wright PetscInt field = 0; // Default field has the coordinates 13427b8c5038SJames Wright 13437b8c5038SJames Wright PetscFunctionBeginUser; 13447b8c5038SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 13457b8c5038SJames Wright PetscCall(DMGetLocalSection(dm_coord, §ion_coord)); 13467b8c5038SJames Wright PetscCall(PetscSectionGetFieldComponents(section_coord, field, num_comp)); 13477b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 13487b8c5038SJames Wright } 13497b8c5038SJames Wright 13507b8c5038SJames Wright /** 13517b8c5038SJames Wright @brief Get the number of components for `dm`'s field 13527b8c5038SJames Wright 13537b8c5038SJames Wright @param[in] dm DM 13547b8c5038SJames Wright @param[in] field The field ID to get the number of components for 13557b8c5038SJames Wright @param[out] num_comp Number of components for the coordinate field 13567b8c5038SJames Wright **/ 13577b8c5038SJames Wright PetscErrorCode DMGetFieldNumComps(DM dm, PetscInt field, PetscInt *num_comp) { 13587b8c5038SJames Wright PetscSection section; 13597b8c5038SJames Wright 13607b8c5038SJames Wright PetscFunctionBeginUser; 13617b8c5038SJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 13627b8c5038SJames Wright PetscCall(PetscSectionGetFieldComponents(section, field, num_comp)); 13637b8c5038SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 13647b8c5038SJames Wright } 1365