Actual source code: ex104.c

  1: static char help[] = "Tests DMPlexCreateColoring().\n\n";

  3: #include <petscdmplex.h>

  5: typedef struct {
  6:   PetscInt depth;
  7:   PetscInt distance;
  8: } AppCtx;

 10: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 11: {
 12:   PetscFunctionBegin;
 13:   options->depth    = 0;
 14:   options->distance = 1;
 15:   PetscOptionsBegin(comm, "", "DMPlexCreateColoring() Test Options", "DMPLEX");
 16:   PetscCall(PetscOptionsInt("-depth", "Stratum depth defining the nodes in the connectivity graph", "ex104.c", options->depth, &options->depth, NULL));
 17:   PetscCall(PetscOptionsInt("-distance", "Coloring distance", "ex104.c", options->distance, &options->distance, NULL));
 18:   PetscOptionsEnd();
 19:   PetscFunctionReturn(PETSC_SUCCESS);
 20: }

 22: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
 23: {
 24:   DM       pdm     = NULL;
 25:   PetscInt overlap = user->distance;
 26:   PetscInt dim;

 28:   PetscFunctionBegin;
 29:   PetscCall(DMCreate(comm, dm));
 30:   PetscCall(DMSetType(*dm, DMPLEX));
 31:   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_TRUE));
 32:   PetscCall(DMSetFromOptions(*dm));
 33:   PetscCall(DMGetDimension(*dm, &dim));
 34:   if (user->depth == dim) {
 35:     PetscCall(DMSetBasicAdjacency(*dm, PETSC_TRUE, PETSC_FALSE));
 36:   } else {
 37:     PetscCall(DMSetBasicAdjacency(*dm, PETSC_FALSE, PETSC_TRUE));
 38:   }
 39:   {
 40:     PetscPartitioner part;
 41:     PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_"));
 42:     PetscCall(DMPlexGetPartitioner(*dm, &part));
 43:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "lb_"));
 44:     PetscCall(PetscPartitionerSetFromOptions(part));
 45:   }
 46:   PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm));
 47:   if (pdm) {
 48:     PetscCall(DMDestroy(dm));
 49:     *dm = pdm;
 50:   }
 51:   PetscFunctionReturn(PETSC_SUCCESS);
 52: }

 54: int main(int argc, char **argv)
 55: {
 56:   DM         dm;
 57:   AppCtx     user;
 58:   PetscInt   ncolors  = 0;
 59:   IS        *iscolors = NULL;
 60:   ISColoring coloring = NULL;

 62:   PetscFunctionBeginUser;
 63:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
 64:   /* Create a BoxMesh */
 65:   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
 66:   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
 67:   /* Color the DMPlex */
 68:   PetscCall(DMPlexCreateColoring(dm, user.depth, user.distance, &coloring));
 69:   PetscCall(ISColoringGetIS(coloring, PETSC_USE_POINTER, &ncolors, &iscolors));
 70:   for (PetscInt c = 0; c < ncolors; c++) {
 71:     PetscCall(ISViewFromOptions(iscolors[c], NULL, "-iscoloring_view"));
 72:   }
 73:   PetscCall(ISColoringRestoreIS(coloring, PETSC_USE_POINTER, &iscolors));
 74:   PetscCall(ISColoringDestroy(&coloring));
 75:   PetscCall(DMDestroy(&dm));
 76:   PetscCall(PetscFinalize());
 77:   return 0;
 78: }

 80: /*TEST

 82:   test:
 83:     nsize: {{1 2}separate output}
 84:     args: -depth {{0 1 2}separate output} -distance {{1 2}separate output} -iscoloring_view -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type simple

 86: TEST*/