12d4ee042Sprj- static char help[] = "Test scalable partitioning on distributed meshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_OVERLAP}; 6c4762a1bSJed Brown 7c4762a1bSJed Brown typedef struct { 8c4762a1bSJed Brown PetscLogEvent createMeshEvent; 9c4762a1bSJed Brown PetscLogStage stages[4]; 10c4762a1bSJed Brown /* Domain and mesh definition */ 11c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 12c4762a1bSJed Brown } AppCtx; 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscErrorCode ierr; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscFunctionBegin; 19c4762a1bSJed Brown options->overlap = PETSC_FALSE; 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 222d4ee042Sprj- ierr = PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex29.c", options->overlap, &options->overlap, NULL,0);CHKERRQ(ierr); 23*1e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 24c4762a1bSJed Brown 25c4762a1bSJed Brown ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshOverlap", &options->stages[STAGE_OVERLAP]);CHKERRQ(ierr); 30c4762a1bSJed Brown PetscFunctionReturn(0); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 34c4762a1bSJed Brown { 35c4762a1bSJed Brown PetscMPIInt rank, size; 36c4762a1bSJed Brown PetscErrorCode ierr; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBegin; 39c4762a1bSJed Brown ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 40ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 41ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 42c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_LOAD]);CHKERRQ(ierr); 4330602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 4430602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 4530602db0SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 47c4762a1bSJed Brown { 48c4762a1bSJed Brown DM pdm = NULL; 49c4762a1bSJed Brown PetscPartitioner part; 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 53c4762a1bSJed Brown /* Distribute mesh over processes */ 54c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 56c4762a1bSJed Brown if (pdm) { 57c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 58c4762a1bSJed Brown *dm = pdm; 59c4762a1bSJed Brown } 60c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr); 6330602db0SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_");CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 6530602db0SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "");CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 67c4762a1bSJed Brown if (user->overlap) { 68c4762a1bSJed Brown DM odm = NULL; 69c4762a1bSJed Brown /* Add the level-1 overlap to refined mesh */ 70c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_OVERLAP]);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = DMPlexDistributeOverlap(*dm, 1, NULL, &odm);CHKERRQ(ierr); 72c4762a1bSJed Brown if (odm) { 73c4762a1bSJed Brown ierr = DMView(odm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 75c4762a1bSJed Brown *dm = odm; 76c4762a1bSJed Brown } 77c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 81c4762a1bSJed Brown PetscFunctionReturn(0); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown int main(int argc, char **argv) 85c4762a1bSJed Brown { 86c4762a1bSJed Brown DM dm, pdm; 87c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 88c4762a1bSJed Brown PetscErrorCode ierr; 89c4762a1bSJed Brown PetscPartitioner part; 90c4762a1bSJed Brown 91c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 92c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = DMPlexDistribute(dm, user.overlap, NULL, &pdm);CHKERRQ(ierr); 97c4762a1bSJed Brown if (pdm) {ierr = DMViewFromOptions(pdm, NULL, "-pdm_view");CHKERRQ(ierr);} 98c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = DMDestroy(&pdm);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = PetscFinalize(); 101c4762a1bSJed Brown return ierr; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /*TEST 105c4762a1bSJed Brown 106c4762a1bSJed Brown test: 107c4762a1bSJed Brown suffix: 0 108c4762a1bSJed Brown requires: ctetgen 10930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -post_dm_refine 2 -petscpartitioner_type simple -dm_view 110c4762a1bSJed Brown test: 111c4762a1bSJed Brown suffix: 1 11230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view 113c4762a1bSJed Brown test: 114c4762a1bSJed Brown suffix: quad_0 115c4762a1bSJed Brown nsize: 2 11630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view 117c4762a1bSJed Brown 118c4762a1bSJed Brown TEST*/ 119