xref: /libCEED/examples/solids/src/setup-dm.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33d8e8822SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53d8e8822SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d8e8822SJeremy L Thompson 
85754ecacSJeremy L Thompson /// @file
95754ecacSJeremy L Thompson /// DM setup for solid mechanics example using PETSc
105754ecacSJeremy L Thompson 
115754ecacSJeremy L Thompson #include "../include/setup-dm.h"
125754ecacSJeremy L Thompson 
13*2b730f8bSJeremy L Thompson #include "../include/boundary.h"
14*2b730f8bSJeremy L Thompson 
155754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
165754ecacSJeremy L Thompson // Setup DM
175754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
185754ecacSJeremy L Thompson PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
195754ecacSJeremy L Thompson   DMLabel label;
205754ecacSJeremy L Thompson 
215754ecacSJeremy L Thompson   PetscFunctionBeginUser;
225754ecacSJeremy L Thompson 
23*2b730f8bSJeremy L Thompson   PetscCall(DMCreateLabel(dm, name));
24*2b730f8bSJeremy L Thompson   PetscCall(DMGetLabel(dm, name, &label));
25*2b730f8bSJeremy L Thompson   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
265754ecacSJeremy L Thompson 
275754ecacSJeremy L Thompson   PetscFunctionReturn(0);
285754ecacSJeremy L Thompson };
295754ecacSJeremy L Thompson 
305754ecacSJeremy L Thompson // Read mesh and distribute DM in parallel
315754ecacSJeremy L Thompson PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx app_ctx, DM *dm) {
325754ecacSJeremy L Thompson   const char      *filename         = app_ctx->mesh_file;
335754ecacSJeremy L Thompson   PetscBool        interpolate      = PETSC_TRUE;
345754ecacSJeremy L Thompson   DM               distributed_mesh = NULL;
355754ecacSJeremy L Thompson   PetscPartitioner part;
365754ecacSJeremy L Thompson 
375754ecacSJeremy L Thompson   PetscFunctionBeginUser;
385754ecacSJeremy L Thompson 
395754ecacSJeremy L Thompson   if (!*filename) {
405754ecacSJeremy L Thompson     PetscInt dim = 3, faces[3] = {3, 3, 3};
41*2b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &dim, NULL));
4220e46440SLeila Ghaffari     if (!dim) dim = 3;
43*2b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, NULL, NULL, NULL, interpolate, dm));
445754ecacSJeremy L Thompson   } else {
45*2b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateFromFile(comm, filename, NULL, interpolate, dm));
465754ecacSJeremy L Thompson   }
475754ecacSJeremy L Thompson 
485754ecacSJeremy L Thompson   // Distribute DM in parallel
49*2b730f8bSJeremy L Thompson   PetscCall(DMPlexGetPartitioner(*dm, &part));
50*2b730f8bSJeremy L Thompson   PetscCall(PetscPartitionerSetFromOptions(part));
51*2b730f8bSJeremy L Thompson   PetscCall(DMPlexDistribute(*dm, 0, NULL, &distributed_mesh));
525754ecacSJeremy L Thompson   if (distributed_mesh) {
53*2b730f8bSJeremy L Thompson     PetscCall(DMDestroy(dm));
545754ecacSJeremy L Thompson     *dm = distributed_mesh;
555754ecacSJeremy L Thompson   }
56*2b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
575754ecacSJeremy L Thompson 
585754ecacSJeremy L Thompson   PetscFunctionReturn(0);
595754ecacSJeremy L Thompson };
605754ecacSJeremy L Thompson 
615754ecacSJeremy L Thompson // Setup DM with FE space of appropriate degree
62*2b730f8bSJeremy L Thompson PetscErrorCode SetupDMByDegree(DM dm, AppCtx app_ctx, PetscInt order, PetscBool boundary, PetscInt num_comp_u) {
6365233696SJeremy L Thompson   MPI_Comm        comm;
645754ecacSJeremy L Thompson   PetscInt        dim;
655754ecacSJeremy L Thompson   PetscFE         fe;
665754ecacSJeremy L Thompson   IS              face_set_is;         // Index Set for Face Sets
675754ecacSJeremy L Thompson   const char     *name = "Face Sets";  // PETSc internal requirement
685754ecacSJeremy L Thompson   PetscInt        num_face_sets;       // Number of FaceSets in face_set_is
695754ecacSJeremy L Thompson   const PetscInt *face_set_ids;        // id of each FaceSet
705754ecacSJeremy L Thompson 
715754ecacSJeremy L Thompson   PetscFunctionBeginUser;
725754ecacSJeremy L Thompson 
735754ecacSJeremy L Thompson   // Setup DM
74*2b730f8bSJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
75*2b730f8bSJeremy L Thompson   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
76*2b730f8bSJeremy L Thompson   PetscCall(PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, order, order, &fe));
77*2b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(dm));
78*2b730f8bSJeremy L Thompson   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
79*2b730f8bSJeremy L Thompson   PetscCall(DMCreateDS(dm));
807ed3e4cdSJeremy L Thompson   {
817ed3e4cdSJeremy L Thompson     /* create FE field for coordinates */
827ed3e4cdSJeremy L Thompson     PetscFE  fe_coords;
837ed3e4cdSJeremy L Thompson     PetscInt num_comp_coord;
84*2b730f8bSJeremy L Thompson     PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
85*2b730f8bSJeremy L Thompson     PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, PETSC_FALSE, 1, 1, &fe_coords));
86*2b730f8bSJeremy L Thompson     PetscCall(DMProjectCoordinates(dm, fe_coords));
87*2b730f8bSJeremy L Thompson     PetscCall(PetscFEDestroy(&fe_coords));
887ed3e4cdSJeremy L Thompson   }
895754ecacSJeremy L Thompson 
905754ecacSJeremy L Thompson   // Add Dirichlet (Essential) boundary
915754ecacSJeremy L Thompson   if (boundary) {
925754ecacSJeremy L Thompson     if (app_ctx->forcing_choice == FORCE_MMS) {
935754ecacSJeremy L Thompson       if (app_ctx->test_mode) {
945754ecacSJeremy L Thompson         // -- Test mode - box mesh
955754ecacSJeremy L Thompson         PetscBool has_label;
965754ecacSJeremy L Thompson         PetscInt  marker_ids[1] = {1};
97*2b730f8bSJeremy L Thompson         PetscCall(DMHasLabel(dm, "marker", &has_label));
985754ecacSJeremy L Thompson         if (!has_label) {
99*2b730f8bSJeremy L Thompson           PetscCall(CreateBCLabel(dm, "marker"));
1005754ecacSJeremy L Thompson         }
1015754ecacSJeremy L Thompson         DMLabel label;
102*2b730f8bSJeremy L Thompson         PetscCall(DMGetLabel(dm, "marker", &label));
103*2b730f8bSJeremy L Thompson         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 1, marker_ids, 0, 0, NULL, (void (*)(void))BCMMS, NULL, NULL, NULL));
1045754ecacSJeremy L Thompson       } else {
1055754ecacSJeremy L Thompson         // -- ExodusII mesh with MMS
106*2b730f8bSJeremy L Thompson         PetscCall(DMGetLabelIdIS(dm, name, &face_set_is));
107*2b730f8bSJeremy L Thompson         PetscCall(ISGetSize(face_set_is, &num_face_sets));
108*2b730f8bSJeremy L Thompson         PetscCall(ISGetIndices(face_set_is, &face_set_ids));
1095754ecacSJeremy L Thompson         DMLabel label;
110*2b730f8bSJeremy L Thompson         PetscCall(DMGetLabel(dm, "Face Sets", &label));
111*2b730f8bSJeremy L Thompson         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, num_face_sets, face_set_ids, 0, 0, NULL, (void (*)(void))BCMMS, NULL, NULL, NULL));
112*2b730f8bSJeremy L Thompson         PetscCall(ISRestoreIndices(face_set_is, &face_set_ids));
113*2b730f8bSJeremy L Thompson         PetscCall(ISDestroy(&face_set_is));
1145754ecacSJeremy L Thompson       }
1155754ecacSJeremy L Thompson     } else {
1167ed3e4cdSJeremy L Thompson       // -- Mesh with user specified BCs
1175754ecacSJeremy L Thompson       DMLabel label;
118*2b730f8bSJeremy L Thompson       PetscCall(DMGetLabel(dm, "Face Sets", &label));
1197ed3e4cdSJeremy L Thompson       // -- Clamp BCs
1205754ecacSJeremy L Thompson       for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) {
1215754ecacSJeremy L Thompson         char bcName[25];
122*2b730f8bSJeremy L Thompson         snprintf(bcName, sizeof bcName, "clamp_%" PetscInt_FMT, app_ctx->bc_clamp_faces[i]);
123*2b730f8bSJeremy L Thompson         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, bcName, label, 1, &app_ctx->bc_clamp_faces[i], 0, 0, NULL, (void (*)(void))BCClamp, NULL,
124*2b730f8bSJeremy L Thompson                                 (void *)&app_ctx->bc_clamp_max[i], NULL));
1255754ecacSJeremy L Thompson       }
1265754ecacSJeremy L Thompson     }
1275754ecacSJeremy L Thompson   }
128*2b730f8bSJeremy L Thompson   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
1295754ecacSJeremy L Thompson 
1305754ecacSJeremy L Thompson   // Cleanup
131*2b730f8bSJeremy L Thompson   PetscCall(PetscFEDestroy(&fe));
1325754ecacSJeremy L Thompson 
1335754ecacSJeremy L Thompson   PetscFunctionReturn(0);
1345754ecacSJeremy L Thompson };
135