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