15754ecacSJeremy L Thompson /// @file 25754ecacSJeremy L Thompson /// DM setup for solid mechanics example using PETSc 35754ecacSJeremy L Thompson 45754ecacSJeremy L Thompson #include "../include/boundary.h" 55754ecacSJeremy L Thompson #include "../include/setup-dm.h" 65754ecacSJeremy L Thompson #include "../include/petsc-macros.h" 75754ecacSJeremy L Thompson 85754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 95754ecacSJeremy L Thompson // Setup DM 105754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 115754ecacSJeremy L Thompson PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 125754ecacSJeremy L Thompson PetscErrorCode ierr; 135754ecacSJeremy L Thompson DMLabel label; 145754ecacSJeremy L Thompson 155754ecacSJeremy L Thompson PetscFunctionBeginUser; 165754ecacSJeremy L Thompson 175754ecacSJeremy L Thompson ierr = DMCreateLabel(dm, name); CHKERRQ(ierr); 185754ecacSJeremy L Thompson ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr); 195754ecacSJeremy L Thompson ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr); 205754ecacSJeremy L Thompson ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr); 215754ecacSJeremy L Thompson 225754ecacSJeremy L Thompson PetscFunctionReturn(0); 235754ecacSJeremy L Thompson }; 245754ecacSJeremy L Thompson 255754ecacSJeremy L Thompson // Read mesh and distribute DM in parallel 265754ecacSJeremy L Thompson PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx app_ctx, DM *dm) { 275754ecacSJeremy L Thompson PetscErrorCode ierr; 285754ecacSJeremy L Thompson const char *filename = app_ctx->mesh_file; 295754ecacSJeremy L Thompson PetscBool interpolate = PETSC_TRUE; 305754ecacSJeremy L Thompson DM distributed_mesh = NULL; 315754ecacSJeremy L Thompson PetscPartitioner part; 325754ecacSJeremy L Thompson 335754ecacSJeremy L Thompson PetscFunctionBeginUser; 345754ecacSJeremy L Thompson 355754ecacSJeremy L Thompson if (!*filename) { 365754ecacSJeremy L Thompson PetscInt dim = 3, faces[3] = {3, 3, 3}; 375754ecacSJeremy L Thompson ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", 385754ecacSJeremy L Thompson faces, &dim, NULL); CHKERRQ(ierr); 39*20e46440SLeila Ghaffari if (!dim) dim = 3; 405754ecacSJeremy L Thompson ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, NULL, 415754ecacSJeremy L Thompson NULL, NULL, interpolate, dm); CHKERRQ(ierr); 425754ecacSJeremy L Thompson } else { 435754ecacSJeremy L Thompson ierr = DMPlexCreateFromFile(comm, filename, interpolate, dm); CHKERRQ(ierr); 445754ecacSJeremy L Thompson } 455754ecacSJeremy L Thompson 465754ecacSJeremy L Thompson // Distribute DM in parallel 475754ecacSJeremy L Thompson ierr = DMPlexGetPartitioner(*dm, &part); CHKERRQ(ierr); 485754ecacSJeremy L Thompson ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 495754ecacSJeremy L Thompson ierr = DMPlexDistribute(*dm, 0, NULL, &distributed_mesh); CHKERRQ(ierr); 505754ecacSJeremy L Thompson if (distributed_mesh) { 515754ecacSJeremy L Thompson ierr = DMDestroy(dm); CHKERRQ(ierr); 525754ecacSJeremy L Thompson *dm = distributed_mesh; 535754ecacSJeremy L Thompson } 545754ecacSJeremy L Thompson ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 555754ecacSJeremy L Thompson 565754ecacSJeremy L Thompson PetscFunctionReturn(0); 575754ecacSJeremy L Thompson }; 585754ecacSJeremy L Thompson 595754ecacSJeremy L Thompson // Create FE by degree 605754ecacSJeremy L Thompson PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc, 615754ecacSJeremy L Thompson PetscBool isSimplex, const char prefix[], 625754ecacSJeremy L Thompson PetscInt order, PetscFE *fem) { 635754ecacSJeremy L Thompson PetscQuadrature q, fq; 645754ecacSJeremy L Thompson DM K; 655754ecacSJeremy L Thompson PetscSpace P; 665754ecacSJeremy L Thompson PetscDualSpace Q; 675754ecacSJeremy L Thompson PetscInt quadPointsPerEdge; 685754ecacSJeremy L Thompson PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 695754ecacSJeremy L Thompson PetscErrorCode ierr; 705754ecacSJeremy L Thompson 715754ecacSJeremy L Thompson PetscFunctionBeginUser; 725754ecacSJeremy L Thompson /* Create space */ 735754ecacSJeremy L Thompson ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P); CHKERRQ(ierr); 745754ecacSJeremy L Thompson ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix); CHKERRQ(ierr); 755754ecacSJeremy L Thompson ierr = PetscSpacePolynomialSetTensor(P, tensor); CHKERRQ(ierr); 765754ecacSJeremy L Thompson ierr = PetscSpaceSetFromOptions(P); CHKERRQ(ierr); 775754ecacSJeremy L Thompson ierr = PetscSpaceSetNumComponents(P, Nc); CHKERRQ(ierr); 785754ecacSJeremy L Thompson ierr = PetscSpaceSetNumVariables(P, dim); CHKERRQ(ierr); 795754ecacSJeremy L Thompson ierr = PetscSpaceSetDegree(P, order, order); CHKERRQ(ierr); 805754ecacSJeremy L Thompson ierr = PetscSpaceSetUp(P); CHKERRQ(ierr); 815754ecacSJeremy L Thompson ierr = PetscSpacePolynomialGetTensor(P, &tensor); CHKERRQ(ierr); 825754ecacSJeremy L Thompson /* Create dual space */ 835754ecacSJeremy L Thompson ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q); 845754ecacSJeremy L Thompson CHKERRQ(ierr); 855754ecacSJeremy L Thompson ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE); CHKERRQ(ierr); 865754ecacSJeremy L Thompson ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix); CHKERRQ(ierr); 875754ecacSJeremy L Thompson ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K); CHKERRQ(ierr); 885754ecacSJeremy L Thompson ierr = PetscDualSpaceSetDM(Q, K); CHKERRQ(ierr); 895754ecacSJeremy L Thompson ierr = DMDestroy(&K); CHKERRQ(ierr); 905754ecacSJeremy L Thompson ierr = PetscDualSpaceSetNumComponents(Q, Nc); CHKERRQ(ierr); 915754ecacSJeremy L Thompson ierr = PetscDualSpaceSetOrder(Q, order); CHKERRQ(ierr); 925754ecacSJeremy L Thompson ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor); CHKERRQ(ierr); 935754ecacSJeremy L Thompson ierr = PetscDualSpaceSetFromOptions(Q); CHKERRQ(ierr); 945754ecacSJeremy L Thompson ierr = PetscDualSpaceSetUp(Q); CHKERRQ(ierr); 955754ecacSJeremy L Thompson /* Create element */ 965754ecacSJeremy L Thompson ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem); CHKERRQ(ierr); 975754ecacSJeremy L Thompson ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix); CHKERRQ(ierr); 985754ecacSJeremy L Thompson ierr = PetscFESetFromOptions(*fem); CHKERRQ(ierr); 995754ecacSJeremy L Thompson ierr = PetscFESetBasisSpace(*fem, P); CHKERRQ(ierr); 1005754ecacSJeremy L Thompson ierr = PetscFESetDualSpace(*fem, Q); CHKERRQ(ierr); 1015754ecacSJeremy L Thompson ierr = PetscFESetNumComponents(*fem, Nc); CHKERRQ(ierr); 1025754ecacSJeremy L Thompson ierr = PetscFESetUp(*fem); CHKERRQ(ierr); 1035754ecacSJeremy L Thompson ierr = PetscSpaceDestroy(&P); CHKERRQ(ierr); 1045754ecacSJeremy L Thompson ierr = PetscDualSpaceDestroy(&Q); CHKERRQ(ierr); 1055754ecacSJeremy L Thompson /* Create quadrature */ 1065754ecacSJeremy L Thompson quadPointsPerEdge = PetscMax(order + 1,1); 1075754ecacSJeremy L Thompson if (isSimplex) { 1085754ecacSJeremy L Thompson ierr = PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, 1095754ecacSJeremy L Thompson &q); CHKERRQ(ierr); 1105754ecacSJeremy L Thompson ierr = PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, 1115754ecacSJeremy L Thompson &fq); CHKERRQ(ierr); 1125754ecacSJeremy L Thompson } else { 1135754ecacSJeremy L Thompson ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, 1145754ecacSJeremy L Thompson &q); CHKERRQ(ierr); 1155754ecacSJeremy L Thompson ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, 1165754ecacSJeremy L Thompson &fq); CHKERRQ(ierr); 1175754ecacSJeremy L Thompson } 1185754ecacSJeremy L Thompson ierr = PetscFESetQuadrature(*fem, q); CHKERRQ(ierr); 1195754ecacSJeremy L Thompson ierr = PetscFESetFaceQuadrature(*fem, fq); CHKERRQ(ierr); 1205754ecacSJeremy L Thompson ierr = PetscQuadratureDestroy(&q); CHKERRQ(ierr); 1215754ecacSJeremy L Thompson ierr = PetscQuadratureDestroy(&fq); CHKERRQ(ierr); 1225754ecacSJeremy L Thompson 1235754ecacSJeremy L Thompson PetscFunctionReturn(0); 1245754ecacSJeremy L Thompson } 1255754ecacSJeremy L Thompson 1265754ecacSJeremy L Thompson // Setup DM with FE space of appropriate degree 1275754ecacSJeremy L Thompson PetscErrorCode SetupDMByDegree(DM dm, AppCtx app_ctx, PetscInt order, 1285754ecacSJeremy L Thompson PetscBool boundary, PetscInt num_comp_u) { 1295754ecacSJeremy L Thompson PetscErrorCode ierr; 1305754ecacSJeremy L Thompson PetscInt dim; 1315754ecacSJeremy L Thompson PetscFE fe; 1325754ecacSJeremy L Thompson IS face_set_is; // Index Set for Face Sets 1335754ecacSJeremy L Thompson const char *name = "Face Sets"; // PETSc internal requirement 1345754ecacSJeremy L Thompson PetscInt num_face_sets; // Number of FaceSets in face_set_is 1355754ecacSJeremy L Thompson const PetscInt *face_set_ids; // id of each FaceSet 1365754ecacSJeremy L Thompson 1375754ecacSJeremy L Thompson PetscFunctionBeginUser; 1385754ecacSJeremy L Thompson 1395754ecacSJeremy L Thompson // Setup DM 1405754ecacSJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 1415754ecacSJeremy L Thompson ierr = PetscFECreateByDegree(dm, dim, num_comp_u, PETSC_FALSE, NULL, order, 1425754ecacSJeremy L Thompson &fe); 1435754ecacSJeremy L Thompson CHKERRQ(ierr); 1445754ecacSJeremy L Thompson ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 1455754ecacSJeremy L Thompson ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 1465754ecacSJeremy L Thompson ierr = DMCreateDS(dm); CHKERRQ(ierr); 1475754ecacSJeremy L Thompson 1485754ecacSJeremy L Thompson // Add Dirichlet (Essential) boundary 1495754ecacSJeremy L Thompson if (boundary) { 1505754ecacSJeremy L Thompson if (app_ctx->forcing_choice == FORCE_MMS) { 1515754ecacSJeremy L Thompson if (app_ctx->test_mode) { 1525754ecacSJeremy L Thompson // -- Test mode - box mesh 1535754ecacSJeremy L Thompson PetscBool has_label; 1545754ecacSJeremy L Thompson PetscInt marker_ids[1] = {1}; 1555754ecacSJeremy L Thompson ierr = DMHasLabel(dm, "marker", &has_label); CHKERRQ(ierr); 1565754ecacSJeremy L Thompson if (!has_label) { 1575754ecacSJeremy L Thompson ierr = CreateBCLabel(dm, "marker"); CHKERRQ(ierr); 1585754ecacSJeremy L Thompson } 1595754ecacSJeremy L Thompson DMLabel label; 1605754ecacSJeremy L Thompson ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr); 1615754ecacSJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, "marker", 1, marker_ids, 1625754ecacSJeremy L Thompson 0, 0, NULL, (void(*)(void))BCMMS, NULL, NULL, NULL); 1635754ecacSJeremy L Thompson CHKERRQ(ierr); 1645754ecacSJeremy L Thompson } else { 1655754ecacSJeremy L Thompson // -- ExodusII mesh with MMS 1665754ecacSJeremy L Thompson ierr = DMGetLabelIdIS(dm, name, &face_set_is); CHKERRQ(ierr); 1675754ecacSJeremy L Thompson ierr = ISGetSize(face_set_is,&num_face_sets); CHKERRQ(ierr); 1685754ecacSJeremy L Thompson ierr = ISGetIndices(face_set_is, &face_set_ids); CHKERRQ(ierr); 1695754ecacSJeremy L Thompson DMLabel label; 1705754ecacSJeremy L Thompson ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 1715754ecacSJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, "Face Sets", 1725754ecacSJeremy L Thompson num_face_sets, face_set_ids, 0, 0, NULL, 1735754ecacSJeremy L Thompson (void(*)(void))BCMMS, NULL, NULL, NULL); 1745754ecacSJeremy L Thompson CHKERRQ(ierr); 1755754ecacSJeremy L Thompson ierr = ISRestoreIndices(face_set_is, &face_set_ids); CHKERRQ(ierr); 1765754ecacSJeremy L Thompson ierr = ISDestroy(&face_set_is); CHKERRQ(ierr); 1775754ecacSJeremy L Thompson } 1785754ecacSJeremy L Thompson } else { 1795754ecacSJeremy L Thompson // -- ExodusII mesh with user specified BCs 1805754ecacSJeremy L Thompson // -- Clamp BCs 1815754ecacSJeremy L Thompson DMLabel label; 1825754ecacSJeremy L Thompson ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 1835754ecacSJeremy L Thompson for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) { 1845754ecacSJeremy L Thompson char bcName[25]; 1855754ecacSJeremy L Thompson snprintf(bcName, sizeof bcName, "clamp_%d", app_ctx->bc_clamp_faces[i]); 1865754ecacSJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, bcName, label, "Face Sets", 1, 1875754ecacSJeremy L Thompson &app_ctx->bc_clamp_faces[i], 0, 0, 1885754ecacSJeremy L Thompson NULL, (void(*)(void))BCClamp, NULL, 1895754ecacSJeremy L Thompson (void *)&app_ctx->bc_clamp_max[i], NULL); 1905754ecacSJeremy L Thompson CHKERRQ(ierr); 1915754ecacSJeremy L Thompson } 1925754ecacSJeremy L Thompson } 1935754ecacSJeremy L Thompson } 1945754ecacSJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 1955754ecacSJeremy L Thompson CHKERRQ(ierr); 1965754ecacSJeremy L Thompson 1975754ecacSJeremy L Thompson // Cleanup 1985754ecacSJeremy L Thompson ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 1995754ecacSJeremy L Thompson 2005754ecacSJeremy L Thompson PetscFunctionReturn(0); 2015754ecacSJeremy L Thompson }; 202