1e83e87a5Sjeremylt #include "../include/petscutils.h" 2e83e87a5Sjeremylt 3e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 4e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType 5e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 69b072555Sjeremylt CeedMemType MemTypeP2C(PetscMemType mem_type) { 79b072555Sjeremylt return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; 8e83e87a5Sjeremylt } 9e83e87a5Sjeremylt 10e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 112c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation 122c58efb6Sjeremylt // ----------------------------------------------------------------------------- 132c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1. Optionally 142c58efb6Sjeremylt // smooth -- see the commented versions at the end. 152c58efb6Sjeremylt static double step(const double a, const double b, double x) { 162c58efb6Sjeremylt if (x <= 0) return a; 172c58efb6Sjeremylt if (x >= 1) return b; 182c58efb6Sjeremylt return a + (b-a) * (x); 192c58efb6Sjeremylt } 202c58efb6Sjeremylt 212c58efb6Sjeremylt // 1D transformation at the right boundary 222c58efb6Sjeremylt static double right(const double eps, const double x) { 232c58efb6Sjeremylt return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1); 242c58efb6Sjeremylt } 252c58efb6Sjeremylt 262c58efb6Sjeremylt // 1D transformation at the left boundary 272c58efb6Sjeremylt static double left(const double eps, const double x) { 282c58efb6Sjeremylt return 1-right(eps,1-x); 292c58efb6Sjeremylt } 302c58efb6Sjeremylt 312c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation 322c58efb6Sjeremylt // The eps parameters are in (0, 1] 332c58efb6Sjeremylt // Uniform mesh is recovered for eps=1 349b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) { 352c58efb6Sjeremylt PetscErrorCode ierr; 362c58efb6Sjeremylt Vec coord; 372c58efb6Sjeremylt PetscInt ncoord; 382c58efb6Sjeremylt PetscScalar *c; 392c58efb6Sjeremylt 402c58efb6Sjeremylt PetscFunctionBeginUser; 419b072555Sjeremylt ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr); 422c58efb6Sjeremylt ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr); 432c58efb6Sjeremylt ierr = VecGetArray(coord, &c); CHKERRQ(ierr); 442c58efb6Sjeremylt 452c58efb6Sjeremylt for (PetscInt i = 0; i < ncoord; i += 3) { 462c58efb6Sjeremylt PetscScalar x = c[i], y = c[i+1], z = c[i+2]; 472c58efb6Sjeremylt PetscInt layer = x*6; 482c58efb6Sjeremylt PetscScalar lambda = (x-layer/6.0)*6; 492c58efb6Sjeremylt c[i] = x; 502c58efb6Sjeremylt 512c58efb6Sjeremylt switch (layer) { 522c58efb6Sjeremylt case 0: 532c58efb6Sjeremylt c[i+1] = left(eps, y); 542c58efb6Sjeremylt c[i+2] = left(eps, z); 552c58efb6Sjeremylt break; 562c58efb6Sjeremylt case 1: 572c58efb6Sjeremylt case 4: 582c58efb6Sjeremylt c[i+1] = step(left(eps, y), right(eps, y), lambda); 592c58efb6Sjeremylt c[i+2] = step(left(eps, z), right(eps, z), lambda); 602c58efb6Sjeremylt break; 612c58efb6Sjeremylt case 2: 622c58efb6Sjeremylt c[i+1] = step(right(eps, y), left(eps, y), lambda/2); 632c58efb6Sjeremylt c[i+2] = step(right(eps, z), left(eps, z), lambda/2); 642c58efb6Sjeremylt break; 652c58efb6Sjeremylt case 3: 662c58efb6Sjeremylt c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2); 672c58efb6Sjeremylt c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2); 682c58efb6Sjeremylt break; 692c58efb6Sjeremylt default: 702c58efb6Sjeremylt c[i+1] = right(eps, y); 712c58efb6Sjeremylt c[i+2] = right(eps, z); 722c58efb6Sjeremylt } 732c58efb6Sjeremylt } 742c58efb6Sjeremylt ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr); 752c58efb6Sjeremylt PetscFunctionReturn(0); 762c58efb6Sjeremylt } 772c58efb6Sjeremylt 782c58efb6Sjeremylt // ----------------------------------------------------------------------------- 79e83e87a5Sjeremylt // Create BC label 80e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 81e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 82751eb813Srezgarshakeri 83e83e87a5Sjeremylt DMLabel label; 84e83e87a5Sjeremylt 85e83e87a5Sjeremylt PetscFunctionBeginUser; 86e83e87a5Sjeremylt 87751eb813Srezgarshakeri PetscCall(DMCreateLabel(dm, name)); 88751eb813Srezgarshakeri PetscCall(DMGetLabel(dm, name, &label)); 89751eb813Srezgarshakeri PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, label)); 90751eb813Srezgarshakeri PetscCall(DMPlexLabelComplete(dm, label)); 91e83e87a5Sjeremylt 92e83e87a5Sjeremylt PetscFunctionReturn(0); 93e83e87a5Sjeremylt }; 94e83e87a5Sjeremylt 95e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 96e83e87a5Sjeremylt // This function sets up a DM for a given degree 97e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 98de1229c5Srezgarshakeri PetscErrorCode SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra, 99*ab1ff315Srezgarshakeri PetscInt num_comp_u, PetscInt dim, bool enforce_bc) { 100e83e87a5Sjeremylt PetscInt ierr, marker_ids[1] = {1}; 101de1229c5Srezgarshakeri PetscInt q_degree = p_degree + q_extra; 102e83e87a5Sjeremylt PetscFE fe; 10365233696SJeremy L Thompson MPI_Comm comm; 104129d8736Srezgarshakeri PetscBool is_simplex = PETSC_TRUE; 105e83e87a5Sjeremylt 106e83e87a5Sjeremylt PetscFunctionBeginUser; 107e83e87a5Sjeremylt 108129d8736Srezgarshakeri // Check if simplex or tensor-product mesh 109129d8736Srezgarshakeri ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr); 110e83e87a5Sjeremylt // Setup FE 11165233696SJeremy L Thompson ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr); 112de1229c5Srezgarshakeri ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree, 113de1229c5Srezgarshakeri q_degree, &fe); CHKERRQ(ierr); 114e83e87a5Sjeremylt ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 115129d8736Srezgarshakeri ierr = DMCreateDS(dm); CHKERRQ(ierr); 116129d8736Srezgarshakeri 117129d8736Srezgarshakeri { 118129d8736Srezgarshakeri // create FE field for coordinates 1197ed3e4cdSJeremy L Thompson PetscFE fe_coords; 1207ed3e4cdSJeremy L Thompson PetscInt num_comp_coord; 1217ed3e4cdSJeremy L Thompson ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 122de1229c5Srezgarshakeri ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree, 1237ed3e4cdSJeremy L Thompson &fe_coords); CHKERRQ(ierr); 1247ed3e4cdSJeremy L Thompson ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 1257ed3e4cdSJeremy L Thompson ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 1267ed3e4cdSJeremy L Thompson } 127de1229c5Srezgarshakeri 128*ab1ff315Srezgarshakeri // Setup Dirichlet BC 129*ab1ff315Srezgarshakeri // Note bp1, bp2 are projection and we don't need to apply BC 130*ab1ff315Srezgarshakeri // For bp3,bp4, the target function is zero on the boundaries 131*ab1ff315Srezgarshakeri // So we pass bcFunc = NULL in DMAddBoundary function 1329b072555Sjeremylt if (enforce_bc) { 1339b072555Sjeremylt PetscBool has_label; 1349b072555Sjeremylt DMHasLabel(dm, "marker", &has_label); 1359b072555Sjeremylt if (!has_label) {CreateBCLabel(dm, "marker");} 13669f5adf1Sjeremylt DMLabel label; 13769f5adf1Sjeremylt ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr); 138b8962995SJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, 139*ab1ff315Srezgarshakeri marker_ids, 0, 0, NULL, NULL, 140b8962995SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 141751eb813Srezgarshakeri PetscCall(DMSetOptionsPrefix(dm, "final_")); 142751eb813Srezgarshakeri PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 143e83e87a5Sjeremylt } 144129d8736Srezgarshakeri 145129d8736Srezgarshakeri if (!is_simplex) { 146129d8736Srezgarshakeri DM dm_coord; 147129d8736Srezgarshakeri ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 148de1229c5Srezgarshakeri ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 149de1229c5Srezgarshakeri CHKERRQ(ierr); 150de1229c5Srezgarshakeri ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); 151de1229c5Srezgarshakeri CHKERRQ(ierr); 152129d8736Srezgarshakeri } 153e83e87a5Sjeremylt ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 154e83e87a5Sjeremylt 155e83e87a5Sjeremylt PetscFunctionReturn(0); 156e83e87a5Sjeremylt }; 157e83e87a5Sjeremylt 158e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 159e83e87a5Sjeremylt // Get CEED restriction data from DMPlex 160e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 1617ed3e4cdSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, 1627ed3e4cdSJeremy L Thompson DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) { 1637ed3e4cdSJeremy L Thompson PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets; 164e83e87a5Sjeremylt PetscErrorCode ierr; 165e83e87a5Sjeremylt 166e83e87a5Sjeremylt PetscFunctionBeginUser; 167e83e87a5Sjeremylt 1687ed3e4cdSJeremy L Thompson ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, 1697ed3e4cdSJeremy L Thompson &elem_size, &num_comp, &num_dof, &elem_restr_offsets); 1707ed3e4cdSJeremy L Thompson CHKERRQ(ierr); 171e83e87a5Sjeremylt 1727ed3e4cdSJeremy L Thompson CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1737ed3e4cdSJeremy L Thompson 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, 1749b072555Sjeremylt elem_restr_offsets, elem_restr); 1759b072555Sjeremylt ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr); 1767ed3e4cdSJeremy L Thompson 177e83e87a5Sjeremylt PetscFunctionReturn(0); 178e83e87a5Sjeremylt }; 179e83e87a5Sjeremylt 180e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 181129d8736Srezgarshakeri // Utility function - convert from DMPolytopeType to CeedElemTopology 182129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 183de1229c5Srezgarshakeri CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) { 184129d8736Srezgarshakeri switch (cell_type) { 185129d8736Srezgarshakeri case DM_POLYTOPE_TRIANGLE: return CEED_TOPOLOGY_TRIANGLE; 186129d8736Srezgarshakeri case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD; 187129d8736Srezgarshakeri case DM_POLYTOPE_TETRAHEDRON: return CEED_TOPOLOGY_TET; 188129d8736Srezgarshakeri case DM_POLYTOPE_HEXAHEDRON: return CEED_TOPOLOGY_HEX; 189129d8736Srezgarshakeri default: return 0; 190129d8736Srezgarshakeri } 191129d8736Srezgarshakeri } 192129d8736Srezgarshakeri 193129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 194f755c37aSrezgarshakeri // Convert DM field to DS field 195129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 196f755c37aSrezgarshakeri PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, 197f755c37aSrezgarshakeri PetscInt *ds_field) { 198129d8736Srezgarshakeri PetscDS ds; 199129d8736Srezgarshakeri IS field_is; 200129d8736Srezgarshakeri const PetscInt *fields; 201129d8736Srezgarshakeri PetscInt num_fields; 202129d8736Srezgarshakeri 203f755c37aSrezgarshakeri PetscFunctionBeginUser; 204f755c37aSrezgarshakeri 205129d8736Srezgarshakeri // Translate dm_field to ds_field 206f755c37aSrezgarshakeri PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds)); 207f755c37aSrezgarshakeri PetscCall(ISGetIndices(field_is, &fields)); 208f755c37aSrezgarshakeri PetscCall(ISGetSize(field_is, &num_fields)); 209129d8736Srezgarshakeri for (PetscInt i = 0; i < num_fields; i++) { 210129d8736Srezgarshakeri if (dm_field == fields[i]) { 211f755c37aSrezgarshakeri *ds_field = i; 212129d8736Srezgarshakeri break; 213129d8736Srezgarshakeri } 214129d8736Srezgarshakeri } 215f755c37aSrezgarshakeri PetscCall(ISRestoreIndices(field_is, &fields)); 216f755c37aSrezgarshakeri 217f755c37aSrezgarshakeri if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, 21851ad7d5bSrezgarshakeri "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 219f755c37aSrezgarshakeri 220f755c37aSrezgarshakeri PetscFunctionReturn(0); 221129d8736Srezgarshakeri } 222129d8736Srezgarshakeri 223f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 224f755c37aSrezgarshakeri // Create libCEED Basis from PetscTabulation 225f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 226f755c37aSrezgarshakeri PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, 227f755c37aSrezgarshakeri PetscInt label_value, PetscInt height, PetscInt face, 228f755c37aSrezgarshakeri PetscFE fe, PetscTabulation basis_tabulation, PetscQuadrature quadrature, 229f755c37aSrezgarshakeri CeedBasis *basis) { 230129d8736Srezgarshakeri 231f755c37aSrezgarshakeri PetscInt first_point; 232129d8736Srezgarshakeri PetscInt ids[1] = {label_value}; 233129d8736Srezgarshakeri DMLabel depth_label; 234129d8736Srezgarshakeri DMPolytopeType cell_type; 235129d8736Srezgarshakeri CeedElemTopology elem_topo; 236f755c37aSrezgarshakeri PetscScalar *q_points, *interp, *grad; 237f755c37aSrezgarshakeri const PetscScalar *q_weights; 238f755c37aSrezgarshakeri PetscDualSpace dual_space; 239f755c37aSrezgarshakeri PetscInt num_dual_basis_vectors; 240f755c37aSrezgarshakeri PetscInt dim, num_comp, P, Q; 241f755c37aSrezgarshakeri 242f755c37aSrezgarshakeri PetscFunctionBeginUser; 243f755c37aSrezgarshakeri 244f755c37aSrezgarshakeri // General basis information 245f755c37aSrezgarshakeri PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 246f755c37aSrezgarshakeri PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 247f755c37aSrezgarshakeri PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 248f755c37aSrezgarshakeri PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 249f755c37aSrezgarshakeri P = num_dual_basis_vectors / num_comp; 250129d8736Srezgarshakeri 251129d8736Srezgarshakeri // Use depth label if no domain label present 252129d8736Srezgarshakeri if (!domain_label) { 253129d8736Srezgarshakeri PetscInt depth; 254129d8736Srezgarshakeri 255f755c37aSrezgarshakeri PetscCall(DMPlexGetDepth(dm, &depth)); 256f755c37aSrezgarshakeri PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 257129d8736Srezgarshakeri ids[0] = depth - height; 258129d8736Srezgarshakeri } 259f755c37aSrezgarshakeri 260129d8736Srezgarshakeri // Get cell interp, grad, and quadrature data 261f755c37aSrezgarshakeri PetscCall(DMGetFirstLabeledPoint(dm, dm, 262f755c37aSrezgarshakeri domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL)); 263f755c37aSrezgarshakeri PetscCall(DMPlexGetCellType(dm, first_point, &cell_type)); 264129d8736Srezgarshakeri elem_topo = ElemTopologyP2C(cell_type); 265129d8736Srezgarshakeri if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, 266129d8736Srezgarshakeri "DMPlex topology not supported"); 267f755c37aSrezgarshakeri { 268f755c37aSrezgarshakeri size_t q_points_size; 269f755c37aSrezgarshakeri const PetscScalar *q_points_petsc; 270f755c37aSrezgarshakeri PetscInt q_dim; 271f755c37aSrezgarshakeri 272f755c37aSrezgarshakeri PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, 273f755c37aSrezgarshakeri &q_weights)); 274f755c37aSrezgarshakeri q_points_size = Q * dim * sizeof(CeedScalar); 275f755c37aSrezgarshakeri PetscCall(PetscCalloc(q_points_size, &q_points)); 276f755c37aSrezgarshakeri for (PetscInt q = 0; q < Q; q++) { 277f755c37aSrezgarshakeri for (PetscInt d = 0; d < q_dim; 278f755c37aSrezgarshakeri d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d]; 279f755c37aSrezgarshakeri } 280f755c37aSrezgarshakeri } 281f755c37aSrezgarshakeri 282129d8736Srezgarshakeri // Convert to libCEED orientation 283f755c37aSrezgarshakeri { 284f755c37aSrezgarshakeri PetscBool is_simplex = PETSC_FALSE; 285f755c37aSrezgarshakeri IS permutation = NULL; 286f755c37aSrezgarshakeri const PetscInt *permutation_indices; 287f755c37aSrezgarshakeri 288f755c37aSrezgarshakeri PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 289f755c37aSrezgarshakeri if (!is_simplex) { 290f755c37aSrezgarshakeri PetscSection section; 291f755c37aSrezgarshakeri 292f755c37aSrezgarshakeri // -- Get permutation 293f755c37aSrezgarshakeri PetscCall(DMGetLocalSection(dm, §ion)); 294f755c37aSrezgarshakeri PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, 295f755c37aSrezgarshakeri num_comp * P, &permutation)); 296f755c37aSrezgarshakeri PetscCall(ISGetIndices(permutation, &permutation_indices)); 297f755c37aSrezgarshakeri } 298f755c37aSrezgarshakeri 299f755c37aSrezgarshakeri // -- Copy interp, grad matrices 300f755c37aSrezgarshakeri PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp)); 301f755c37aSrezgarshakeri PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad)); 302129d8736Srezgarshakeri const CeedInt c = 0; 303129d8736Srezgarshakeri for (CeedInt q = 0; q < Q; q++) { 304f755c37aSrezgarshakeri for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) { 305f755c37aSrezgarshakeri CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed 306f755c37aSrezgarshakeri * num_comp]; 307f755c37aSrezgarshakeri 308f755c37aSrezgarshakeri interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + 309f755c37aSrezgarshakeri p_petsc) * num_comp + c]; 310129d8736Srezgarshakeri for (CeedInt d = 0; d < dim; d++) { 311f755c37aSrezgarshakeri grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][((( 312f755c37aSrezgarshakeri face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d]; 313129d8736Srezgarshakeri } 314129d8736Srezgarshakeri } 315129d8736Srezgarshakeri } 316f755c37aSrezgarshakeri 317f755c37aSrezgarshakeri // -- Cleanup 318f755c37aSrezgarshakeri if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices)); 319f755c37aSrezgarshakeri PetscCall(ISDestroy(&permutation)); 320f755c37aSrezgarshakeri } 321f755c37aSrezgarshakeri 322f755c37aSrezgarshakeri // Finally, create libCEED basis 323f755c37aSrezgarshakeri CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, 324f755c37aSrezgarshakeri q_weights, basis); 325f755c37aSrezgarshakeri PetscCall(PetscFree(q_points)); 326f755c37aSrezgarshakeri PetscCall(PetscFree(interp)); 327f755c37aSrezgarshakeri PetscCall(PetscFree(grad)); 328f755c37aSrezgarshakeri 329f755c37aSrezgarshakeri PetscFunctionReturn(0); 330f755c37aSrezgarshakeri } 331f755c37aSrezgarshakeri 332f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 333f755c37aSrezgarshakeri // Get CEED Basis from DMPlex 334f755c37aSrezgarshakeri // ----------------------------------------------------------------------------- 335f755c37aSrezgarshakeri PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, 336f755c37aSrezgarshakeri CeedInt label_value, CeedInt height, 337d4d45553Srezgarshakeri CeedInt dm_field, BPData bp_data, CeedBasis *basis) { 338f755c37aSrezgarshakeri PetscDS ds; 339f755c37aSrezgarshakeri PetscFE fe; 340f755c37aSrezgarshakeri PetscQuadrature quadrature; 341f755c37aSrezgarshakeri PetscBool is_simplex = PETSC_TRUE; 342f755c37aSrezgarshakeri PetscInt ds_field = -1; 343f755c37aSrezgarshakeri 344f755c37aSrezgarshakeri PetscFunctionBeginUser; 345f755c37aSrezgarshakeri 346f755c37aSrezgarshakeri // Get element information 347f755c37aSrezgarshakeri PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds)); 348f755c37aSrezgarshakeri PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 349f755c37aSrezgarshakeri PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 350f755c37aSrezgarshakeri PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 351f755c37aSrezgarshakeri PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 352f755c37aSrezgarshakeri 353f755c37aSrezgarshakeri // Check if simplex or tensor-product mesh 354f755c37aSrezgarshakeri PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 355f755c37aSrezgarshakeri 356f755c37aSrezgarshakeri // Build libCEED basis 357f755c37aSrezgarshakeri if (is_simplex) { 358f755c37aSrezgarshakeri PetscTabulation basis_tabulation; 359f755c37aSrezgarshakeri PetscInt num_derivatives = 1, face = 0; 360f755c37aSrezgarshakeri 361f755c37aSrezgarshakeri PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation)); 362f755c37aSrezgarshakeri PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, 363f755c37aSrezgarshakeri face, fe, basis_tabulation, quadrature, basis)); 364129d8736Srezgarshakeri } else { 365f755c37aSrezgarshakeri PetscDualSpace dual_space; 366f755c37aSrezgarshakeri PetscInt num_dual_basis_vectors; 367f755c37aSrezgarshakeri PetscInt dim, num_comp, P, Q; 368f755c37aSrezgarshakeri 369f755c37aSrezgarshakeri PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 370f755c37aSrezgarshakeri PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 371f755c37aSrezgarshakeri PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 372f755c37aSrezgarshakeri PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 373f755c37aSrezgarshakeri P = num_dual_basis_vectors / num_comp; 374f755c37aSrezgarshakeri PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL)); 375f755c37aSrezgarshakeri 376129d8736Srezgarshakeri CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim)); 377129d8736Srezgarshakeri CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim)); 378129d8736Srezgarshakeri 379d4d45553Srezgarshakeri CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, 380d4d45553Srezgarshakeri bp_data.q_mode, basis); 381129d8736Srezgarshakeri } 382129d8736Srezgarshakeri 383129d8736Srezgarshakeri PetscFunctionReturn(0); 384f755c37aSrezgarshakeri } 385129d8736Srezgarshakeri 386129d8736Srezgarshakeri // ----------------------------------------------------------------------------- 387de1229c5Srezgarshakeri // Utilities 388de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 389de1229c5Srezgarshakeri 390de1229c5Srezgarshakeri // Utility function, compute three factors of an integer 391de1229c5Srezgarshakeri static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 392de1229c5Srezgarshakeri for (PetscInt d=0, size_left=size; d<3; d++) { 393de1229c5Srezgarshakeri PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d))); 394de1229c5Srezgarshakeri while (try * (size_left / try) != size_left) try++; 395de1229c5Srezgarshakeri m[reverse ? 2-d : d] = try; 396de1229c5Srezgarshakeri size_left /= try; 397de1229c5Srezgarshakeri } 398de1229c5Srezgarshakeri } 399de1229c5Srezgarshakeri 400de1229c5Srezgarshakeri static int Max3(const PetscInt a[3]) { 401de1229c5Srezgarshakeri return PetscMax(a[0], PetscMax(a[1], a[2])); 402de1229c5Srezgarshakeri } 403de1229c5Srezgarshakeri 404de1229c5Srezgarshakeri static int Min3(const PetscInt a[3]) { 405de1229c5Srezgarshakeri return PetscMin(a[0], PetscMin(a[1], a[2])); 406de1229c5Srezgarshakeri } 407de1229c5Srezgarshakeri 408de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 409de1229c5Srezgarshakeri // Create distribute dm 410de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 411de1229c5Srezgarshakeri PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) { 412de1229c5Srezgarshakeri PetscErrorCode ierr; 413de1229c5Srezgarshakeri 414de1229c5Srezgarshakeri PetscFunctionBeginUser; 415de1229c5Srezgarshakeri // Setup DM 416de1229c5Srezgarshakeri if (rp->read_mesh) { 417de1229c5Srezgarshakeri ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE, 418de1229c5Srezgarshakeri dm); 419de1229c5Srezgarshakeri CHKERRQ(ierr); 420de1229c5Srezgarshakeri } else { 421de1229c5Srezgarshakeri if (rp->user_l_nodes) { 422de1229c5Srezgarshakeri // Find a nicely composite number of elements no less than global nodes 423de1229c5Srezgarshakeri PetscMPIInt size; 424de1229c5Srezgarshakeri ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr); 425de1229c5Srezgarshakeri for (PetscInt g_elem = 426de1229c5Srezgarshakeri PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim)); 427de1229c5Srezgarshakeri ; 428de1229c5Srezgarshakeri g_elem++) { 429de1229c5Srezgarshakeri Split3(g_elem, rp->mesh_elem, true); 430de1229c5Srezgarshakeri if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break; 431de1229c5Srezgarshakeri } 432de1229c5Srezgarshakeri } 433*ab1ff315Srezgarshakeri 434de1229c5Srezgarshakeri ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex, 435de1229c5Srezgarshakeri rp->mesh_elem, 436de1229c5Srezgarshakeri NULL, NULL, NULL, PETSC_TRUE, dm); CHKERRQ(ierr); 437de1229c5Srezgarshakeri } 438de1229c5Srezgarshakeri 439de1229c5Srezgarshakeri ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 440de1229c5Srezgarshakeri ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 441de1229c5Srezgarshakeri 442de1229c5Srezgarshakeri PetscFunctionReturn(0); 443de1229c5Srezgarshakeri } 444de1229c5Srezgarshakeri 445de1229c5Srezgarshakeri // ----------------------------------------------------------------------------- 446