12d4ee042Sprj- static char help[] = "Test scalable partitioning on distributed meshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 59371c9d4SSatish Balay enum { 69371c9d4SSatish Balay STAGE_LOAD, 79371c9d4SSatish Balay STAGE_DISTRIBUTE, 89371c9d4SSatish Balay STAGE_REFINE, 99371c9d4SSatish Balay STAGE_OVERLAP 109371c9d4SSatish Balay }; 11c4762a1bSJed Brown 12c4762a1bSJed Brown typedef struct { 13c4762a1bSJed Brown PetscLogEvent createMeshEvent; 14c4762a1bSJed Brown PetscLogStage stages[4]; 15c4762a1bSJed Brown /* Domain and mesh definition */ 16c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 17c4762a1bSJed Brown } AppCtx; 18c4762a1bSJed Brown 19d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 20d71ae5a4SJacob Faibussowitsch { 21c4762a1bSJed Brown PetscFunctionBegin; 22c4762a1bSJed Brown options->overlap = PETSC_FALSE; 23c4762a1bSJed Brown 24d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex29.c", options->overlap, &options->overlap, NULL, 0)); 26d0609cedSBarry Smith PetscOptionsEnd(); 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent)); 299566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD])); 309566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE])); 319566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE])); 329566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshOverlap", &options->stages[STAGE_OVERLAP])); 33*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 37d71ae5a4SJacob Faibussowitsch { 38c4762a1bSJed Brown PetscMPIInt rank, size; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscFunctionBegin; 419566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(user->createMeshEvent, 0, 0, 0, 0)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 449566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD])); 459566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 469566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 479566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 489566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 49c4762a1bSJed Brown { 50c4762a1bSJed Brown DM pdm = NULL; 51c4762a1bSJed Brown PetscPartitioner part; 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 549566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 55c4762a1bSJed Brown /* Distribute mesh over processes */ 569566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE])); 579566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm)); 58c4762a1bSJed Brown if (pdm) { 599566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 60c4762a1bSJed Brown *dm = pdm; 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 63c4762a1bSJed Brown } 649566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE])); 659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "post_")); 669566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "")); 689566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 69c4762a1bSJed Brown if (user->overlap) { 70c4762a1bSJed Brown DM odm = NULL; 71c4762a1bSJed Brown /* Add the level-1 overlap to refined mesh */ 729566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_OVERLAP])); 739566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeOverlap(*dm, 1, NULL, &odm)); 74c4762a1bSJed Brown if (odm) { 759566063dSJacob Faibussowitsch PetscCall(DMView(odm, PETSC_VIEWER_STDOUT_WORLD)); 769566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 77c4762a1bSJed Brown *dm = odm; 78c4762a1bSJed Brown } 799566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 80c4762a1bSJed Brown } 819566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(user->createMeshEvent, 0, 0, 0, 0)); 83*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 87d71ae5a4SJacob Faibussowitsch { 88c4762a1bSJed Brown DM dm, pdm; 89c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 90c4762a1bSJed Brown PetscPartitioner part; 91c4762a1bSJed Brown 92327415f7SBarry Smith PetscFunctionBeginUser; 939566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 949566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 959566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 969566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 979566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 989566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, user.overlap, NULL, &pdm)); 999566063dSJacob Faibussowitsch if (pdm) PetscCall(DMViewFromOptions(pdm, NULL, "-pdm_view")); 1009566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1019566063dSJacob Faibussowitsch PetscCall(DMDestroy(&pdm)); 1029566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 103b122ec5aSJacob Faibussowitsch return 0; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /*TEST 107c4762a1bSJed Brown 108c4762a1bSJed Brown test: 109c4762a1bSJed Brown suffix: 0 110c4762a1bSJed Brown requires: ctetgen 11130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -post_dm_refine 2 -petscpartitioner_type simple -dm_view 112c4762a1bSJed Brown test: 113c4762a1bSJed Brown suffix: 1 11430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view 115c4762a1bSJed Brown test: 116c4762a1bSJed Brown suffix: quad_0 117c4762a1bSJed Brown nsize: 2 11830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view 119c4762a1bSJed Brown 120c4762a1bSJed Brown TEST*/ 121