1 /// @file 2 /// DM setup for solid mechanics example using PETSc 3 4 #include "../include/boundary.h" 5 #include "../include/setup-dm.h" 6 7 // ----------------------------------------------------------------------------- 8 // Setup DM 9 // ----------------------------------------------------------------------------- 10 PetscErrorCode CreateBCLabel(DM dm, const char name[]) { 11 PetscErrorCode ierr; 12 DMLabel label; 13 14 PetscFunctionBeginUser; 15 16 ierr = DMCreateLabel(dm, name); CHKERRQ(ierr); 17 ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr); 18 ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr); 19 ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr); 20 21 PetscFunctionReturn(0); 22 }; 23 24 // Read mesh and distribute DM in parallel 25 PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx app_ctx, DM *dm) { 26 PetscErrorCode ierr; 27 const char *filename = app_ctx->mesh_file; 28 PetscBool interpolate = PETSC_TRUE; 29 DM distributed_mesh = NULL; 30 PetscPartitioner part; 31 32 PetscFunctionBeginUser; 33 34 if (!*filename) { 35 PetscInt dim = 3, faces[3] = {3, 3, 3}; 36 ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", 37 faces, &dim, NULL); CHKERRQ(ierr); 38 if (!dim) dim = 3; 39 ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, NULL, 40 NULL, NULL, interpolate, dm); CHKERRQ(ierr); 41 } else { 42 ierr = DMPlexCreateFromFile(comm, filename, NULL, interpolate, dm); 43 CHKERRQ(ierr); 44 } 45 46 // Distribute DM in parallel 47 ierr = DMPlexGetPartitioner(*dm, &part); CHKERRQ(ierr); 48 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 49 ierr = DMPlexDistribute(*dm, 0, NULL, &distributed_mesh); CHKERRQ(ierr); 50 if (distributed_mesh) { 51 ierr = DMDestroy(dm); CHKERRQ(ierr); 52 *dm = distributed_mesh; 53 } 54 ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 55 56 PetscFunctionReturn(0); 57 }; 58 59 // Setup DM with FE space of appropriate degree 60 PetscErrorCode SetupDMByDegree(DM dm, AppCtx app_ctx, PetscInt order, 61 PetscBool boundary, PetscInt num_comp_u) { 62 PetscErrorCode ierr; 63 MPI_Comm comm; 64 PetscInt dim; 65 PetscFE fe; 66 IS face_set_is; // Index Set for Face Sets 67 const char *name = "Face Sets"; // PETSc internal requirement 68 PetscInt num_face_sets; // Number of FaceSets in face_set_is 69 const PetscInt *face_set_ids; // id of each FaceSet 70 71 PetscFunctionBeginUser; 72 73 // Setup DM 74 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 75 ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr); 76 ierr = PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, order, order, 77 &fe); CHKERRQ(ierr); 78 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 79 ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 80 ierr = DMCreateDS(dm); CHKERRQ(ierr); 81 { 82 /* create FE field for coordinates */ 83 PetscFE fe_coords; 84 PetscInt num_comp_coord; 85 ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 86 ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, PETSC_FALSE, 1, 1, 87 &fe_coords); CHKERRQ(ierr); 88 ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 89 ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 90 } 91 92 // Add Dirichlet (Essential) boundary 93 if (boundary) { 94 if (app_ctx->forcing_choice == FORCE_MMS) { 95 if (app_ctx->test_mode) { 96 // -- Test mode - box mesh 97 PetscBool has_label; 98 PetscInt marker_ids[1] = {1}; 99 ierr = DMHasLabel(dm, "marker", &has_label); CHKERRQ(ierr); 100 if (!has_label) { 101 ierr = CreateBCLabel(dm, "marker"); CHKERRQ(ierr); 102 } 103 DMLabel label; 104 ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr); 105 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 1, marker_ids, 106 0, 0, NULL, (void(*)(void))BCMMS, NULL, NULL, NULL); 107 CHKERRQ(ierr); 108 } else { 109 // -- ExodusII mesh with MMS 110 ierr = DMGetLabelIdIS(dm, name, &face_set_is); CHKERRQ(ierr); 111 ierr = ISGetSize(face_set_is,&num_face_sets); CHKERRQ(ierr); 112 ierr = ISGetIndices(face_set_is, &face_set_ids); CHKERRQ(ierr); 113 DMLabel label; 114 ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 115 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 116 num_face_sets, face_set_ids, 0, 0, NULL, 117 (void(*)(void))BCMMS, NULL, NULL, NULL); 118 CHKERRQ(ierr); 119 ierr = ISRestoreIndices(face_set_is, &face_set_ids); CHKERRQ(ierr); 120 ierr = ISDestroy(&face_set_is); CHKERRQ(ierr); 121 } 122 } else { 123 // -- Mesh with user specified BCs 124 DMLabel label; 125 ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 126 // -- Clamp BCs 127 for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) { 128 char bcName[25]; 129 snprintf(bcName, sizeof bcName, "clamp_%d", app_ctx->bc_clamp_faces[i]); 130 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, bcName, label, 1, 131 &app_ctx->bc_clamp_faces[i], 0, 0, 132 NULL, (void(*)(void))BCClamp, NULL, 133 (void *)&app_ctx->bc_clamp_max[i], NULL); 134 CHKERRQ(ierr); 135 } 136 } 137 } 138 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 139 CHKERRQ(ierr); 140 141 // Cleanup 142 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 143 144 PetscFunctionReturn(0); 145 }; 146