xref: /petsc/src/dm/impls/plex/tests/ex70.c (revision aee58fc643cde3a74e29e886bde44205d937e3a7)
141e3e744Sksagiyam static char help[] = "Tests for DMPlexMarkBoundaryFaces()\n\n";
241e3e744Sksagiyam 
341e3e744Sksagiyam #include <petscdmplex.h>
441e3e744Sksagiyam #include <petscsf.h>
541e3e744Sksagiyam 
641e3e744Sksagiyam typedef struct {
741e3e744Sksagiyam   PetscInt overlap; /* The overlap size used when partitioning */
841e3e744Sksagiyam } AppCtx;
941e3e744Sksagiyam 
ProcessOptions(MPI_Comm comm,AppCtx * options)1041e3e744Sksagiyam PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
1141e3e744Sksagiyam {
1241e3e744Sksagiyam   PetscFunctionBegin;
1341e3e744Sksagiyam   options->overlap = 0;
1441e3e744Sksagiyam 
1541e3e744Sksagiyam   PetscOptionsBegin(comm, "", "Options for DMPlexMarkBoundaryFaces() problem", "DMPLEX");
1641e3e744Sksagiyam   PetscCall(PetscOptionsBoundedInt("-overlap", "The overlap size used", "ex70.c", options->overlap, &options->overlap, NULL, 0));
1741e3e744Sksagiyam   PetscOptionsEnd();
1841e3e744Sksagiyam   PetscFunctionReturn(PETSC_SUCCESS);
1941e3e744Sksagiyam }
2041e3e744Sksagiyam 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)2141e3e744Sksagiyam static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
2241e3e744Sksagiyam {
2341e3e744Sksagiyam   PetscFunctionBeginUser;
2441e3e744Sksagiyam   {
2541e3e744Sksagiyam     const PetscInt faces[2] = {2, 2};
2641e3e744Sksagiyam 
27*42108689Sksagiyam     PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, 0, PETSC_TRUE, dm));
2841e3e744Sksagiyam   }
2941e3e744Sksagiyam   {
3041e3e744Sksagiyam     PetscPartitioner part;
3141e3e744Sksagiyam     PetscInt        *sizes  = NULL;
3241e3e744Sksagiyam     PetscInt        *points = NULL;
3341e3e744Sksagiyam     PetscMPIInt      rank, size;
3441e3e744Sksagiyam 
3541e3e744Sksagiyam     PetscCallMPI(MPI_Comm_rank(comm, &rank));
3641e3e744Sksagiyam     PetscCallMPI(MPI_Comm_size(comm, &size));
3741e3e744Sksagiyam     if (rank == 0) {
3841e3e744Sksagiyam       PetscInt sizes1[2]  = {4, 4};
3941e3e744Sksagiyam       PetscInt points1[8] = {3, 5, 6, 7, 0, 1, 2, 4};
4041e3e744Sksagiyam 
4141e3e744Sksagiyam       PetscCall(PetscMalloc2(2, &sizes, 8, &points));
4241e3e744Sksagiyam       PetscCall(PetscArraycpy(sizes, sizes1, 2));
4341e3e744Sksagiyam       PetscCall(PetscArraycpy(points, points1, 8));
4441e3e744Sksagiyam     }
4541e3e744Sksagiyam     PetscCall(DMPlexGetPartitioner(*dm, &part));
4641e3e744Sksagiyam     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
4741e3e744Sksagiyam     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
4841e3e744Sksagiyam     PetscCall(PetscFree2(sizes, points));
4941e3e744Sksagiyam   }
5041e3e744Sksagiyam   {
5141e3e744Sksagiyam     DM dmDist = NULL;
5241e3e744Sksagiyam 
5341e3e744Sksagiyam     PetscCall(DMSetAdjacency(*dm, -1, PETSC_FALSE, PETSC_TRUE));
5441e3e744Sksagiyam     PetscCall(DMPlexDistribute(*dm, user->overlap, NULL, &dmDist));
5541e3e744Sksagiyam     if (dmDist) {
5641e3e744Sksagiyam       PetscCall(DMDestroy(dm));
5741e3e744Sksagiyam       *dm = dmDist;
5841e3e744Sksagiyam     }
5941e3e744Sksagiyam   }
6041e3e744Sksagiyam   PetscFunctionReturn(PETSC_SUCCESS);
6141e3e744Sksagiyam }
6241e3e744Sksagiyam 
main(int argc,char ** argv)6341e3e744Sksagiyam int main(int argc, char **argv)
6441e3e744Sksagiyam {
6541e3e744Sksagiyam   DM          dm;
6641e3e744Sksagiyam   DMLabel     extLabel;
6741e3e744Sksagiyam   MPI_Comm    comm;
6841e3e744Sksagiyam   PetscMPIInt size;
6941e3e744Sksagiyam   AppCtx      user;
7041e3e744Sksagiyam 
7141e3e744Sksagiyam   PetscFunctionBeginUser;
7241e3e744Sksagiyam   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
7341e3e744Sksagiyam   comm = PETSC_COMM_WORLD;
7441e3e744Sksagiyam   PetscCallMPI(MPI_Comm_size(comm, &size));
7541e3e744Sksagiyam   if (size != 2) {
7641e3e744Sksagiyam     PetscCall(PetscPrintf(comm, "This example is specifically designed for comm size == 2.\n"));
7741e3e744Sksagiyam     PetscCall(PetscFinalize());
7841e3e744Sksagiyam     return 0;
7941e3e744Sksagiyam   }
8041e3e744Sksagiyam   PetscCall(ProcessOptions(comm, &user));
8141e3e744Sksagiyam   PetscCall(CreateMesh(comm, &user, &dm));
8241e3e744Sksagiyam   PetscCall(DMCreateLabel(dm, "exterior_facets"));
8341e3e744Sksagiyam   PetscCall(DMGetLabel(dm, "exterior_facets", &extLabel));
8441e3e744Sksagiyam   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, extLabel));
8541e3e744Sksagiyam   PetscCall(PetscObjectSetName((PetscObject)dm, "Example_DM"));
8641e3e744Sksagiyam   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
8741e3e744Sksagiyam   PetscCall(DMDestroy(&dm));
8841e3e744Sksagiyam   PetscCall(PetscFinalize());
8941e3e744Sksagiyam   return 0;
9041e3e744Sksagiyam }
9141e3e744Sksagiyam 
9241e3e744Sksagiyam /*TEST
9341e3e744Sksagiyam 
9441e3e744Sksagiyam   test:
9541e3e744Sksagiyam     nsize: 2
9641e3e744Sksagiyam     requires: triangle
9741e3e744Sksagiyam     args: -overlap {{0 1}separate output} -dm_view ascii::ascii_info_detail
9841e3e744Sksagiyam 
9941e3e744Sksagiyam TEST*/
100