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