xref: /libCEED/examples/solids/src/setup-dm.c (revision caee03026e6576cbf7a399c2fc51bb918c77f451)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// DM setup for solid mechanics example using PETSc
10 
11 #include "../include/setup-dm.h"
12 
13 #include "../include/boundary.h"
14 
15 // -----------------------------------------------------------------------------
16 // Setup DM
17 // -----------------------------------------------------------------------------
18 PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
19   DMLabel label;
20 
21   PetscFunctionBeginUser;
22 
23   PetscCall(DMCreateLabel(dm, name));
24   PetscCall(DMGetLabel(dm, name, &label));
25   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
26 
27   PetscFunctionReturn(0);
28 };
29 
30 // Read mesh and distribute DM in parallel
31 PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx app_ctx, DM *dm) {
32   const char      *filename         = app_ctx->mesh_file;
33   PetscBool        interpolate      = PETSC_TRUE;
34   DM               distributed_mesh = NULL;
35   PetscPartitioner part;
36 
37   PetscFunctionBeginUser;
38 
39   if (!*filename) {
40     PetscInt dim = 3, faces[3] = {3, 3, 3};
41     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &dim, NULL));
42     if (!dim) dim = 3;
43     PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, NULL, NULL, NULL, interpolate, dm));
44   } else {
45     PetscCall(DMPlexCreateFromFile(comm, filename, NULL, interpolate, dm));
46   }
47 
48   // Distribute DM in parallel
49   PetscCall(DMPlexGetPartitioner(*dm, &part));
50   PetscCall(PetscPartitionerSetFromOptions(part));
51   PetscCall(DMPlexDistribute(*dm, 0, NULL, &distributed_mesh));
52   if (distributed_mesh) {
53     PetscCall(DMDestroy(dm));
54     *dm = distributed_mesh;
55   }
56   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
57 
58   PetscFunctionReturn(0);
59 };
60 
61 // Setup DM with FE space of appropriate degree
62 PetscErrorCode SetupDMByDegree(DM dm, AppCtx app_ctx, PetscInt order, PetscBool boundary, PetscInt num_comp_u) {
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   PetscCall(DMGetDimension(dm, &dim));
75   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
76   PetscCall(PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, order, order, &fe));
77   PetscCall(DMSetFromOptions(dm));
78   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
79   PetscCall(DMCreateDS(dm));
80   {
81     /* create FE field for coordinates */
82     PetscFE  fe_coords;
83     PetscInt num_comp_coord;
84     PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
85     PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, PETSC_FALSE, 1, 1, &fe_coords));
86     PetscCall(DMProjectCoordinates(dm, fe_coords));
87     PetscCall(PetscFEDestroy(&fe_coords));
88   }
89 
90   // Add Dirichlet (Essential) boundary
91   if (boundary) {
92     if (app_ctx->forcing_choice == FORCE_MMS) {
93       if (app_ctx->test_mode) {
94         // -- Test mode - box mesh
95         PetscBool has_label;
96         PetscInt  marker_ids[1] = {1};
97         PetscCall(DMHasLabel(dm, "marker", &has_label));
98         if (!has_label) {
99           PetscCall(CreateBCLabel(dm, "marker"));
100         }
101         DMLabel label;
102         PetscCall(DMGetLabel(dm, "marker", &label));
103         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 1, marker_ids, 0, 0, NULL, (void (*)(void))BCMMS, NULL, NULL, NULL));
104       } else {
105         // -- ExodusII mesh with MMS
106         PetscCall(DMGetLabelIdIS(dm, name, &face_set_is));
107         PetscCall(ISGetSize(face_set_is, &num_face_sets));
108         PetscCall(ISGetIndices(face_set_is, &face_set_ids));
109         DMLabel label;
110         PetscCall(DMGetLabel(dm, "Face Sets", &label));
111         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, num_face_sets, face_set_ids, 0, 0, NULL, (void (*)(void))BCMMS, NULL, NULL, NULL));
112         PetscCall(ISRestoreIndices(face_set_is, &face_set_ids));
113         PetscCall(ISDestroy(&face_set_is));
114       }
115     } else {
116       // -- Mesh with user specified BCs
117       DMLabel label;
118       PetscCall(DMGetLabel(dm, "Face Sets", &label));
119       // -- Clamp BCs
120       for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) {
121         char bcName[25];
122         snprintf(bcName, sizeof bcName, "clamp_%" PetscInt_FMT, app_ctx->bc_clamp_faces[i]);
123         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, bcName, label, 1, &app_ctx->bc_clamp_faces[i], 0, 0, NULL, (void (*)(void))BCClamp, NULL,
124                                 (void *)&app_ctx->bc_clamp_max[i], NULL));
125       }
126     }
127   }
128   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
129 
130   // Cleanup
131   PetscCall(PetscFEDestroy(&fe));
132 
133   PetscFunctionReturn(0);
134 };
135