1e83e87a5Sjeremylt #include "../include/petscutils.h" 2e83e87a5Sjeremylt 3e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 4e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType 5e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 62b730f8bSJeremy L Thompson CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; } 7e83e87a5Sjeremylt 8e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 92c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation 102c58efb6Sjeremylt // ----------------------------------------------------------------------------- 112c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1. Optionally 122c58efb6Sjeremylt // smooth -- see the commented versions at the end. 132c58efb6Sjeremylt static double step(const double a, const double b, double x) { 142c58efb6Sjeremylt if (x <= 0) return a; 152c58efb6Sjeremylt if (x >= 1) return b; 162c58efb6Sjeremylt return a + (b - a) * (x); 172c58efb6Sjeremylt } 182c58efb6Sjeremylt 192c58efb6Sjeremylt // 1D transformation at the right boundary 202b730f8bSJeremy L Thompson static double right(const double eps, const double x) { return (x <= 0.5) ? (2 - eps) * x : 1 + eps * (x - 1); } 212c58efb6Sjeremylt 222c58efb6Sjeremylt // 1D transformation at the left boundary 232b730f8bSJeremy L Thompson static double left(const double eps, const double x) { return 1 - right(eps, 1 - x); } 242c58efb6Sjeremylt 252c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation 262c58efb6Sjeremylt // The eps parameters are in (0, 1] 272c58efb6Sjeremylt // Uniform mesh is recovered for eps=1 289b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) { 292c58efb6Sjeremylt Vec coord; 302c58efb6Sjeremylt PetscInt ncoord; 312c58efb6Sjeremylt PetscScalar *c; 322c58efb6Sjeremylt 332c58efb6Sjeremylt PetscFunctionBeginUser; 342b730f8bSJeremy L Thompson 352b730f8bSJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm_orig, &coord)); 362b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(coord, &ncoord)); 372b730f8bSJeremy L Thompson PetscCall(VecGetArray(coord, &c)); 382c58efb6Sjeremylt 392c58efb6Sjeremylt for (PetscInt i = 0; i < ncoord; i += 3) { 402c58efb6Sjeremylt PetscScalar x = c[i], y = c[i + 1], z = c[i + 2]; 412c58efb6Sjeremylt PetscInt layer = x * 6; 422c58efb6Sjeremylt PetscScalar lambda = (x - layer / 6.0) * 6; 432c58efb6Sjeremylt c[i] = x; 442c58efb6Sjeremylt 452c58efb6Sjeremylt switch (layer) { 462c58efb6Sjeremylt case 0: 472c58efb6Sjeremylt c[i + 1] = left(eps, y); 482c58efb6Sjeremylt c[i + 2] = left(eps, z); 492c58efb6Sjeremylt break; 502c58efb6Sjeremylt case 1: 512c58efb6Sjeremylt case 4: 522c58efb6Sjeremylt c[i + 1] = step(left(eps, y), right(eps, y), lambda); 532c58efb6Sjeremylt c[i + 2] = step(left(eps, z), right(eps, z), lambda); 542c58efb6Sjeremylt break; 552c58efb6Sjeremylt case 2: 562c58efb6Sjeremylt c[i + 1] = step(right(eps, y), left(eps, y), lambda / 2); 572c58efb6Sjeremylt c[i + 2] = step(right(eps, z), left(eps, z), lambda / 2); 582c58efb6Sjeremylt break; 592c58efb6Sjeremylt case 3: 602c58efb6Sjeremylt c[i + 1] = step(right(eps, y), left(eps, y), (1 + lambda) / 2); 612c58efb6Sjeremylt c[i + 2] = step(right(eps, z), left(eps, z), (1 + lambda) / 2); 622c58efb6Sjeremylt break; 632c58efb6Sjeremylt default: 642c58efb6Sjeremylt c[i + 1] = right(eps, y); 652c58efb6Sjeremylt c[i + 2] = right(eps, z); 662c58efb6Sjeremylt } 672c58efb6Sjeremylt } 682b730f8bSJeremy L Thompson PetscCall(VecRestoreArray(coord, &c)); 69ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 702c58efb6Sjeremylt } 712c58efb6Sjeremylt 722c58efb6Sjeremylt // ----------------------------------------------------------------------------- 73e83e87a5Sjeremylt // Create BC label 74e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 75e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 76e83e87a5Sjeremylt DMLabel label; 77e83e87a5Sjeremylt 78e83e87a5Sjeremylt PetscFunctionBeginUser; 79e83e87a5Sjeremylt 80751eb813Srezgarshakeri PetscCall(DMCreateLabel(dm, name)); 81751eb813Srezgarshakeri PetscCall(DMGetLabel(dm, name, &label)); 82751eb813Srezgarshakeri PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, label)); 83751eb813Srezgarshakeri PetscCall(DMPlexLabelComplete(dm, label)); 84e83e87a5Sjeremylt 85ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 86fd8f24faSSebastian Grimberg } 87e83e87a5Sjeremylt 88e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 89e83e87a5Sjeremylt // This function sets up a DM for a given degree 90e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 912b730f8bSJeremy L Thompson PetscErrorCode SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra, PetscInt num_comp_u, PetscInt dim, bool enforce_bc) { 922b730f8bSJeremy L Thompson PetscInt marker_ids[1] = {1}; 93de1229c5Srezgarshakeri PetscInt q_degree = p_degree + q_extra; 94e83e87a5Sjeremylt PetscFE fe; 9565233696SJeremy L Thompson MPI_Comm comm; 96129d8736Srezgarshakeri PetscBool is_simplex = PETSC_TRUE; 97e83e87a5Sjeremylt 98e83e87a5Sjeremylt PetscFunctionBeginUser; 99e83e87a5Sjeremylt 100129d8736Srezgarshakeri // Check if simplex or tensor-product mesh 1012b730f8bSJeremy L Thompson PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 102e83e87a5Sjeremylt // Setup FE 1032b730f8bSJeremy L Thompson PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1042b730f8bSJeremy L Thompson PetscCall(PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree, q_degree, &fe)); 1052b730f8bSJeremy L Thompson PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 1062b730f8bSJeremy L Thompson PetscCall(DMCreateDS(dm)); 107129d8736Srezgarshakeri 108129d8736Srezgarshakeri { 109129d8736Srezgarshakeri // create FE field for coordinates 1107ed3e4cdSJeremy L Thompson PetscFE fe_coords; 1117ed3e4cdSJeremy L Thompson PetscInt num_comp_coord; 1122b730f8bSJeremy L Thompson PetscCall(DMGetCoordinateDim(dm, &num_comp_coord)); 1132b730f8bSJeremy L Thompson PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree, &fe_coords)); 114*49a40c8aSKenneth E. Jansen PetscCall(DMSetCoordinateDisc(dm, fe_coords, PETSC_TRUE)); 1152b730f8bSJeremy L Thompson PetscCall(PetscFEDestroy(&fe_coords)); 1167ed3e4cdSJeremy L Thompson } 117de1229c5Srezgarshakeri 118ab1ff315Srezgarshakeri // Setup Dirichlet BC 119ab1ff315Srezgarshakeri // Note bp1, bp2 are projection and we don't need to apply BC 120ab1ff315Srezgarshakeri // For bp3,bp4, the target function is zero on the boundaries 121ab1ff315Srezgarshakeri // So we pass bcFunc = NULL in DMAddBoundary function 1229b072555Sjeremylt if (enforce_bc) { 1239b072555Sjeremylt PetscBool has_label; 1249b072555Sjeremylt DMHasLabel(dm, "marker", &has_label); 1252b730f8bSJeremy L Thompson if (!has_label) { 1262b730f8bSJeremy L Thompson CreateBCLabel(dm, "marker"); 1272b730f8bSJeremy L Thompson } 12869f5adf1Sjeremylt DMLabel label; 1292b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, "marker", &label)); 1302b730f8bSJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, marker_ids, 0, 0, NULL, NULL, NULL, NULL, NULL)); 131751eb813Srezgarshakeri PetscCall(DMSetOptionsPrefix(dm, "final_")); 132751eb813Srezgarshakeri PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 133e83e87a5Sjeremylt } 134129d8736Srezgarshakeri 135129d8736Srezgarshakeri if (!is_simplex) { 136129d8736Srezgarshakeri DM dm_coord; 1372b730f8bSJeremy L Thompson PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1382b730f8bSJeremy L Thompson PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 1392b730f8bSJeremy L Thompson PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 140129d8736Srezgarshakeri } 1412b730f8bSJeremy L Thompson PetscCall(PetscFEDestroy(&fe)); 142e83e87a5Sjeremylt 143ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 144fd8f24faSSebastian Grimberg } 145e83e87a5Sjeremylt 146e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 147e83e87a5Sjeremylt // Get CEED restriction data from DMPlex 148e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 1492b730f8bSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) { 1507ed3e4cdSJeremy L Thompson PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets; 151e83e87a5Sjeremylt 152e83e87a5Sjeremylt PetscFunctionBeginUser; 153e83e87a5Sjeremylt 1542b730f8bSJeremy L Thompson PetscCall(DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets)); 155e83e87a5Sjeremylt 1562b730f8bSJeremy L Thompson CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets, elem_restr); 1572b730f8bSJeremy L Thompson PetscCall(PetscFree(elem_restr_offsets)); 1587ed3e4cdSJeremy L Thompson 159ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 160fd8f24faSSebastian Grimberg } 161e83e87a5Sjeremylt 162e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 163129d8736Srezgarshakeri // Utility function - convert from DMPolytopeType to CeedElemTopology 164129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 165de1229c5Srezgarshakeri CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) { 166129d8736Srezgarshakeri switch (cell_type) { 1672b730f8bSJeremy L Thompson case DM_POLYTOPE_TRIANGLE: 1682b730f8bSJeremy L Thompson return CEED_TOPOLOGY_TRIANGLE; 1692b730f8bSJeremy L Thompson case DM_POLYTOPE_QUADRILATERAL: 1702b730f8bSJeremy L Thompson return CEED_TOPOLOGY_QUAD; 1712b730f8bSJeremy L Thompson case DM_POLYTOPE_TETRAHEDRON: 1722b730f8bSJeremy L Thompson return CEED_TOPOLOGY_TET; 1732b730f8bSJeremy L Thompson case DM_POLYTOPE_HEXAHEDRON: 1742b730f8bSJeremy L Thompson return CEED_TOPOLOGY_HEX; 1752b730f8bSJeremy L Thompson default: 1762b730f8bSJeremy L Thompson return 0; 177129d8736Srezgarshakeri } 178129d8736Srezgarshakeri } 179129d8736Srezgarshakeri 180129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 181f755c37aSrezgarshakeri // Convert DM field to DS field 182129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 1832b730f8bSJeremy L Thompson PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) { 184129d8736Srezgarshakeri PetscDS ds; 185129d8736Srezgarshakeri IS field_is; 186129d8736Srezgarshakeri const PetscInt *fields; 187129d8736Srezgarshakeri PetscInt num_fields; 188129d8736Srezgarshakeri 189f755c37aSrezgarshakeri PetscFunctionBeginUser; 190f755c37aSrezgarshakeri 191129d8736Srezgarshakeri // Translate dm_field to ds_field 192baf96a30SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL)); 193f755c37aSrezgarshakeri PetscCall(ISGetIndices(field_is, &fields)); 194f755c37aSrezgarshakeri PetscCall(ISGetSize(field_is, &num_fields)); 195129d8736Srezgarshakeri for (PetscInt i = 0; i < num_fields; i++) { 196129d8736Srezgarshakeri if (dm_field == fields[i]) { 197f755c37aSrezgarshakeri *ds_field = i; 198129d8736Srezgarshakeri break; 199129d8736Srezgarshakeri } 200129d8736Srezgarshakeri } 201f755c37aSrezgarshakeri PetscCall(ISRestoreIndices(field_is, &fields)); 202f755c37aSrezgarshakeri 2032b730f8bSJeremy L Thompson if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 204f755c37aSrezgarshakeri 205ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 206129d8736Srezgarshakeri } 207129d8736Srezgarshakeri 208f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 209f755c37aSrezgarshakeri // Create libCEED Basis from PetscTabulation 210f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 2112b730f8bSJeremy L Thompson PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe, 2122b730f8bSJeremy L Thompson PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) { 213f755c37aSrezgarshakeri PetscInt first_point; 214129d8736Srezgarshakeri PetscInt ids[1] = {label_value}; 215129d8736Srezgarshakeri DMLabel depth_label; 216129d8736Srezgarshakeri DMPolytopeType cell_type; 217129d8736Srezgarshakeri CeedElemTopology elem_topo; 218f755c37aSrezgarshakeri PetscScalar *q_points, *interp, *grad; 219f755c37aSrezgarshakeri const PetscScalar *q_weights; 220f755c37aSrezgarshakeri PetscDualSpace dual_space; 221f755c37aSrezgarshakeri PetscInt num_dual_basis_vectors; 222f755c37aSrezgarshakeri PetscInt dim, num_comp, P, Q; 223f755c37aSrezgarshakeri 224f755c37aSrezgarshakeri PetscFunctionBeginUser; 225f755c37aSrezgarshakeri 226f755c37aSrezgarshakeri // General basis information 227f755c37aSrezgarshakeri PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 228f755c37aSrezgarshakeri PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 229f755c37aSrezgarshakeri PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 230f755c37aSrezgarshakeri PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 231f755c37aSrezgarshakeri P = num_dual_basis_vectors / num_comp; 232129d8736Srezgarshakeri 233129d8736Srezgarshakeri // Use depth label if no domain label present 234129d8736Srezgarshakeri if (!domain_label) { 235129d8736Srezgarshakeri PetscInt depth; 236129d8736Srezgarshakeri 237f755c37aSrezgarshakeri PetscCall(DMPlexGetDepth(dm, &depth)); 238f755c37aSrezgarshakeri PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 239129d8736Srezgarshakeri ids[0] = depth - height; 240129d8736Srezgarshakeri } 241f755c37aSrezgarshakeri 242129d8736Srezgarshakeri // Get cell interp, grad, and quadrature data 2432b730f8bSJeremy L Thompson PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL)); 244f755c37aSrezgarshakeri PetscCall(DMPlexGetCellType(dm, first_point, &cell_type)); 245129d8736Srezgarshakeri elem_topo = ElemTopologyP2C(cell_type); 2462b730f8bSJeremy L Thompson if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported"); 247f755c37aSrezgarshakeri { 248f755c37aSrezgarshakeri size_t q_points_size; 249f755c37aSrezgarshakeri const PetscScalar *q_points_petsc; 250f755c37aSrezgarshakeri PetscInt q_dim; 251f755c37aSrezgarshakeri 2522b730f8bSJeremy L Thompson PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights)); 253f755c37aSrezgarshakeri q_points_size = Q * dim * sizeof(CeedScalar); 254f755c37aSrezgarshakeri PetscCall(PetscCalloc(q_points_size, &q_points)); 255f755c37aSrezgarshakeri for (PetscInt q = 0; q < Q; q++) { 2562b730f8bSJeremy L Thompson for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d]; 257f755c37aSrezgarshakeri } 258f755c37aSrezgarshakeri } 259f755c37aSrezgarshakeri 260129d8736Srezgarshakeri // Convert to libCEED orientation 261f755c37aSrezgarshakeri { 262f755c37aSrezgarshakeri PetscBool is_simplex = PETSC_FALSE; 263f755c37aSrezgarshakeri IS permutation = NULL; 264f755c37aSrezgarshakeri const PetscInt *permutation_indices; 265f755c37aSrezgarshakeri 266f755c37aSrezgarshakeri PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 267f755c37aSrezgarshakeri if (!is_simplex) { 268f755c37aSrezgarshakeri PetscSection section; 269f755c37aSrezgarshakeri 270f755c37aSrezgarshakeri // -- Get permutation 271f755c37aSrezgarshakeri PetscCall(DMGetLocalSection(dm, §ion)); 2722b730f8bSJeremy L Thompson PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation)); 273f755c37aSrezgarshakeri PetscCall(ISGetIndices(permutation, &permutation_indices)); 274f755c37aSrezgarshakeri } 275f755c37aSrezgarshakeri 276f755c37aSrezgarshakeri // -- Copy interp, grad matrices 277f755c37aSrezgarshakeri PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp)); 278f755c37aSrezgarshakeri PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad)); 279129d8736Srezgarshakeri const CeedInt c = 0; 280129d8736Srezgarshakeri for (CeedInt q = 0; q < Q; q++) { 281f755c37aSrezgarshakeri for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) { 2822b730f8bSJeremy L Thompson CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp]; 283f755c37aSrezgarshakeri 2842b730f8bSJeremy L Thompson interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c]; 285129d8736Srezgarshakeri for (CeedInt d = 0; d < dim; d++) { 2862b730f8bSJeremy L Thompson grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d]; 287129d8736Srezgarshakeri } 288129d8736Srezgarshakeri } 289129d8736Srezgarshakeri } 290f755c37aSrezgarshakeri 291f755c37aSrezgarshakeri // -- Cleanup 292f755c37aSrezgarshakeri if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices)); 293f755c37aSrezgarshakeri PetscCall(ISDestroy(&permutation)); 294f755c37aSrezgarshakeri } 295f755c37aSrezgarshakeri 296f755c37aSrezgarshakeri // Finally, create libCEED basis 2972b730f8bSJeremy L Thompson CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis); 298f755c37aSrezgarshakeri PetscCall(PetscFree(q_points)); 299f755c37aSrezgarshakeri PetscCall(PetscFree(interp)); 300f755c37aSrezgarshakeri PetscCall(PetscFree(grad)); 301f755c37aSrezgarshakeri 302ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 303f755c37aSrezgarshakeri } 304f755c37aSrezgarshakeri 305f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 306f755c37aSrezgarshakeri // Get CEED Basis from DMPlex 307f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 3082b730f8bSJeremy L Thompson PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, BPData bp_data, 3092b730f8bSJeremy L Thompson CeedBasis *basis) { 310f755c37aSrezgarshakeri PetscDS ds; 311f755c37aSrezgarshakeri PetscFE fe; 312f755c37aSrezgarshakeri PetscQuadrature quadrature; 313f755c37aSrezgarshakeri PetscBool is_simplex = PETSC_TRUE; 314f755c37aSrezgarshakeri PetscInt ds_field = -1; 315f755c37aSrezgarshakeri 316f755c37aSrezgarshakeri PetscFunctionBeginUser; 317f755c37aSrezgarshakeri 318f755c37aSrezgarshakeri // Get element information 319baf96a30SJames Wright PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 320f755c37aSrezgarshakeri PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 321f755c37aSrezgarshakeri PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 322f755c37aSrezgarshakeri PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 323f755c37aSrezgarshakeri PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 324f755c37aSrezgarshakeri 325f755c37aSrezgarshakeri // Check if simplex or tensor-product mesh 326f755c37aSrezgarshakeri PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 327f755c37aSrezgarshakeri 328f755c37aSrezgarshakeri // Build libCEED basis 329f755c37aSrezgarshakeri if (is_simplex) { 330f755c37aSrezgarshakeri PetscTabulation basis_tabulation; 331f755c37aSrezgarshakeri PetscInt num_derivatives = 1, face = 0; 332f755c37aSrezgarshakeri 333f755c37aSrezgarshakeri PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation)); 3342b730f8bSJeremy L Thompson PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis)); 335129d8736Srezgarshakeri } else { 336f755c37aSrezgarshakeri PetscDualSpace dual_space; 337f755c37aSrezgarshakeri PetscInt num_dual_basis_vectors; 338f755c37aSrezgarshakeri PetscInt dim, num_comp, P, Q; 339f755c37aSrezgarshakeri 340f755c37aSrezgarshakeri PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 341f755c37aSrezgarshakeri PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 342f755c37aSrezgarshakeri PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 343f755c37aSrezgarshakeri PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 344f755c37aSrezgarshakeri P = num_dual_basis_vectors / num_comp; 345f755c37aSrezgarshakeri PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL)); 346f755c37aSrezgarshakeri 347129d8736Srezgarshakeri CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim)); 348129d8736Srezgarshakeri CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim)); 349129d8736Srezgarshakeri 3502b730f8bSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, bp_data.q_mode, basis); 351129d8736Srezgarshakeri } 352129d8736Srezgarshakeri 353ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 354f755c37aSrezgarshakeri } 355129d8736Srezgarshakeri 356129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 357de1229c5Srezgarshakeri // Utilities 358de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 359de1229c5Srezgarshakeri 360de1229c5Srezgarshakeri // Utility function, compute three factors of an integer 361de1229c5Srezgarshakeri static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 362de1229c5Srezgarshakeri for (PetscInt d = 0, size_left = size; d < 3; d++) { 363de1229c5Srezgarshakeri PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1. / (3 - d))); 364de1229c5Srezgarshakeri while (try * (size_left / try) != size_left) try++; 365de1229c5Srezgarshakeri m[reverse ? 2 - d : d] = try; 366de1229c5Srezgarshakeri size_left /= try; 367de1229c5Srezgarshakeri } 368de1229c5Srezgarshakeri } 369de1229c5Srezgarshakeri 3702b730f8bSJeremy L Thompson static int Max3(const PetscInt a[3]) { return PetscMax(a[0], PetscMax(a[1], a[2])); } 371de1229c5Srezgarshakeri 3722b730f8bSJeremy L Thompson static int Min3(const PetscInt a[3]) { return PetscMin(a[0], PetscMin(a[1], a[2])); } 373de1229c5Srezgarshakeri 374de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 375de1229c5Srezgarshakeri // Create distribute dm 376de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 377de1229c5Srezgarshakeri PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) { 378de1229c5Srezgarshakeri PetscFunctionBeginUser; 379de1229c5Srezgarshakeri // Setup DM 380de1229c5Srezgarshakeri if (rp->read_mesh) { 3812b730f8bSJeremy L Thompson PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE, dm)); 382de1229c5Srezgarshakeri } else { 383de1229c5Srezgarshakeri if (rp->user_l_nodes) { 384de1229c5Srezgarshakeri // Find a nicely composite number of elements no less than global nodes 385de1229c5Srezgarshakeri PetscMPIInt size; 3862b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(rp->comm, &size)); 3872b730f8bSJeremy L Thompson for (PetscInt g_elem = PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));; g_elem++) { 388de1229c5Srezgarshakeri Split3(g_elem, rp->mesh_elem, true); 389de1229c5Srezgarshakeri if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break; 390de1229c5Srezgarshakeri } 391de1229c5Srezgarshakeri } 392ab1ff315Srezgarshakeri 3932b730f8bSJeremy L Thompson PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex, rp->mesh_elem, NULL, NULL, NULL, PETSC_TRUE, dm)); 394de1229c5Srezgarshakeri } 395de1229c5Srezgarshakeri 3962b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(*dm)); 3972b730f8bSJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 398de1229c5Srezgarshakeri 399ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 400de1229c5Srezgarshakeri } 401de1229c5Srezgarshakeri 402de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 403