1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdmplex.h>
1549aac155SJeremy L Thompson
162b730f8bSJeremy L Thompson #include "../include/boundary.h"
172b730f8bSJeremy L Thompson
185754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
195754ecacSJeremy L Thompson // Setup DM
205754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
CreateBCLabel(DM dm,const char name[])215754ecacSJeremy L Thompson PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
225754ecacSJeremy L Thompson DMLabel label;
235754ecacSJeremy L Thompson
245754ecacSJeremy L Thompson PetscFunctionBeginUser;
255754ecacSJeremy L Thompson
262b730f8bSJeremy L Thompson PetscCall(DMCreateLabel(dm, name));
272b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, name, &label));
282b730f8bSJeremy L Thompson PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
295754ecacSJeremy L Thompson
30ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
315754ecacSJeremy L Thompson };
325754ecacSJeremy L Thompson
335754ecacSJeremy L Thompson // Read mesh and distribute DM in parallel
CreateDistributedDM(MPI_Comm comm,AppCtx app_ctx,DM * dm)345754ecacSJeremy L Thompson PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx app_ctx, DM *dm) {
355754ecacSJeremy L Thompson const char *filename = app_ctx->mesh_file;
365754ecacSJeremy L Thompson PetscBool interpolate = PETSC_TRUE;
375754ecacSJeremy L Thompson DM distributed_mesh = NULL;
385754ecacSJeremy L Thompson PetscPartitioner part;
395754ecacSJeremy L Thompson
405754ecacSJeremy L Thompson PetscFunctionBeginUser;
415754ecacSJeremy L Thompson
425754ecacSJeremy L Thompson if (!*filename) {
435754ecacSJeremy L Thompson PetscInt dim = 3, faces[3] = {3, 3, 3};
442b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &dim, NULL));
4520e46440SLeila Ghaffari if (!dim) dim = 3;
467cf95199SJeremy L Thompson PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, NULL, NULL, NULL, interpolate, 0, PETSC_FALSE, dm));
475754ecacSJeremy L Thompson } else {
482b730f8bSJeremy L Thompson PetscCall(DMPlexCreateFromFile(comm, filename, NULL, interpolate, dm));
495754ecacSJeremy L Thompson }
505754ecacSJeremy L Thompson
515754ecacSJeremy L Thompson // Distribute DM in parallel
522b730f8bSJeremy L Thompson PetscCall(DMPlexGetPartitioner(*dm, &part));
532b730f8bSJeremy L Thompson PetscCall(PetscPartitionerSetFromOptions(part));
542b730f8bSJeremy L Thompson PetscCall(DMPlexDistribute(*dm, 0, NULL, &distributed_mesh));
555754ecacSJeremy L Thompson if (distributed_mesh) {
562b730f8bSJeremy L Thompson PetscCall(DMDestroy(dm));
575754ecacSJeremy L Thompson *dm = distributed_mesh;
585754ecacSJeremy L Thompson }
592b730f8bSJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
605754ecacSJeremy L Thompson
61ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
625754ecacSJeremy L Thompson };
635754ecacSJeremy L Thompson
645754ecacSJeremy L Thompson // Setup DM with FE space of appropriate degree
SetupDMByDegree(DM dm,AppCtx app_ctx,PetscInt order,PetscBool boundary,PetscInt num_comp_u)652b730f8bSJeremy L Thompson PetscErrorCode SetupDMByDegree(DM dm, AppCtx app_ctx, PetscInt order, PetscBool boundary, PetscInt num_comp_u) {
6665233696SJeremy L Thompson MPI_Comm comm;
675754ecacSJeremy L Thompson PetscInt dim;
685754ecacSJeremy L Thompson PetscFE fe;
695754ecacSJeremy L Thompson IS face_set_is; // Index Set for Face Sets
705754ecacSJeremy L Thompson const char *name = "Face Sets"; // PETSc internal requirement
715754ecacSJeremy L Thompson PetscInt num_face_sets; // Number of FaceSets in face_set_is
725754ecacSJeremy L Thompson const PetscInt *face_set_ids; // id of each FaceSet
735754ecacSJeremy L Thompson
745754ecacSJeremy L Thompson PetscFunctionBeginUser;
755754ecacSJeremy L Thompson
765754ecacSJeremy L Thompson // Setup DM
772b730f8bSJeremy L Thompson PetscCall(DMGetDimension(dm, &dim));
782b730f8bSJeremy L Thompson PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
792b730f8bSJeremy L Thompson PetscCall(PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, order, order, &fe));
802b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(dm));
812b730f8bSJeremy L Thompson PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
822b730f8bSJeremy L Thompson PetscCall(DMCreateDS(dm));
837ed3e4cdSJeremy L Thompson {
847ed3e4cdSJeremy L Thompson /* create FE field for coordinates */
857ed3e4cdSJeremy L Thompson PetscFE fe_coords;
867ed3e4cdSJeremy L Thompson PetscInt num_comp_coord;
872b730f8bSJeremy L Thompson PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
882b730f8bSJeremy L Thompson PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, PETSC_FALSE, 1, 1, &fe_coords));
8949a40c8aSKenneth E. Jansen PetscCall(DMSetCoordinateDisc(dm, fe_coords, PETSC_TRUE));
902b730f8bSJeremy L Thompson PetscCall(PetscFEDestroy(&fe_coords));
917ed3e4cdSJeremy L Thompson }
925754ecacSJeremy L Thompson
935754ecacSJeremy L Thompson // Add Dirichlet (Essential) boundary
945754ecacSJeremy L Thompson if (boundary) {
955754ecacSJeremy L Thompson if (app_ctx->forcing_choice == FORCE_MMS) {
965754ecacSJeremy L Thompson if (app_ctx->test_mode) {
975754ecacSJeremy L Thompson // -- Test mode - box mesh
985754ecacSJeremy L Thompson PetscBool has_label;
995754ecacSJeremy L Thompson PetscInt marker_ids[1] = {1};
1002b730f8bSJeremy L Thompson PetscCall(DMHasLabel(dm, "marker", &has_label));
1015754ecacSJeremy L Thompson if (!has_label) {
1022b730f8bSJeremy L Thompson PetscCall(CreateBCLabel(dm, "marker"));
1035754ecacSJeremy L Thompson }
1045754ecacSJeremy L Thompson DMLabel label;
1052b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, "marker", &label));
1062b730f8bSJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 1, marker_ids, 0, 0, NULL, (void (*)(void))BCMMS, NULL, NULL, NULL));
1075754ecacSJeremy L Thompson } else {
1085754ecacSJeremy L Thompson // -- ExodusII mesh with MMS
1092b730f8bSJeremy L Thompson PetscCall(DMGetLabelIdIS(dm, name, &face_set_is));
1102b730f8bSJeremy L Thompson PetscCall(ISGetSize(face_set_is, &num_face_sets));
1112b730f8bSJeremy L Thompson PetscCall(ISGetIndices(face_set_is, &face_set_ids));
1125754ecacSJeremy L Thompson DMLabel label;
1132b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label));
1142b730f8bSJeremy 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));
1152b730f8bSJeremy L Thompson PetscCall(ISRestoreIndices(face_set_is, &face_set_ids));
1162b730f8bSJeremy L Thompson PetscCall(ISDestroy(&face_set_is));
1175754ecacSJeremy L Thompson }
1185754ecacSJeremy L Thompson } else {
1197ed3e4cdSJeremy L Thompson // -- Mesh with user specified BCs
1205754ecacSJeremy L Thompson DMLabel label;
1212b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label));
1227ed3e4cdSJeremy L Thompson // -- Clamp BCs
1235754ecacSJeremy L Thompson for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) {
1245754ecacSJeremy L Thompson char bcName[25];
1252b730f8bSJeremy L Thompson snprintf(bcName, sizeof bcName, "clamp_%" PetscInt_FMT, app_ctx->bc_clamp_faces[i]);
1262b730f8bSJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, bcName, label, 1, &app_ctx->bc_clamp_faces[i], 0, 0, NULL, (void (*)(void))BCClamp, NULL,
1272b730f8bSJeremy L Thompson (void *)&app_ctx->bc_clamp_max[i], NULL));
1285754ecacSJeremy L Thompson }
1295754ecacSJeremy L Thompson }
1305754ecacSJeremy L Thompson }
1312b730f8bSJeremy L Thompson PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
1325754ecacSJeremy L Thompson
1335754ecacSJeremy L Thompson // Cleanup
1342b730f8bSJeremy L Thompson PetscCall(PetscFEDestroy(&fe));
1355754ecacSJeremy L Thompson
136ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
1375754ecacSJeremy L Thompson };
138