1*2d4ee042Sprj- 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 dim; /* The topological mesh dimension */ 12c4762a1bSJed Brown PetscInt faces[3]; /* Number of faces per dimension */ 13c4762a1bSJed Brown PetscBool simplex; /* Use simplices or hexes */ 14c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 15c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 16c4762a1bSJed Brown } AppCtx; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown PetscInt dim; 21c4762a1bSJed Brown PetscErrorCode ierr; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 24c4762a1bSJed Brown options->dim = 2; 25c4762a1bSJed Brown options->simplex = PETSC_TRUE; 26c4762a1bSJed Brown options->filename[0] = '\0'; 27c4762a1bSJed Brown options->overlap = PETSC_FALSE; 28c4762a1bSJed Brown options->faces[0] = 1; 29c4762a1bSJed Brown options->faces[1] = 1; 30c4762a1bSJed Brown options->faces[2] = 1; 31c4762a1bSJed Brown 32c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 33*2d4ee042Sprj- ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex29.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 34*2d4ee042Sprj- ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex29.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = PetscOptionsString("-filename", "The mesh file", "", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr); 36*2d4ee042Sprj- ierr = PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex29.c", options->overlap, &options->overlap, NULL,0);CHKERRQ(ierr); 37c4762a1bSJed Brown dim = options->dim; 38*2d4ee042Sprj- ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex29.c", options->faces, &dim, NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown if (dim) options->dim = dim; 40c4762a1bSJed Brown ierr = PetscOptionsEnd(); 41c4762a1bSJed Brown 42c4762a1bSJed Brown ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshOverlap", &options->stages[STAGE_OVERLAP]);CHKERRQ(ierr); 47c4762a1bSJed Brown PetscFunctionReturn(0); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown PetscInt dim = user->dim; 53c4762a1bSJed Brown PetscInt *faces = user->faces; 54c4762a1bSJed Brown PetscBool simplex = user->simplex; 55c4762a1bSJed Brown const char *filename = user->filename; 56c4762a1bSJed Brown size_t len; 57c4762a1bSJed Brown PetscMPIInt rank, size; 58c4762a1bSJed Brown PetscErrorCode ierr; 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscFunctionBegin; 61c4762a1bSJed Brown ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_LOAD]);CHKERRQ(ierr); 66c4762a1bSJed Brown if (len) { 67c4762a1bSJed Brown ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr); 68c4762a1bSJed Brown } else { 69c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 72c4762a1bSJed Brown { 73c4762a1bSJed Brown DM pdm = NULL; 74c4762a1bSJed Brown PetscPartitioner part; 75c4762a1bSJed Brown 76c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 78c4762a1bSJed Brown /* Distribute mesh over processes */ 79c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 81c4762a1bSJed Brown if (pdm) { 82c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 83c4762a1bSJed Brown *dm = pdm; 84c4762a1bSJed Brown } 85c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 90c4762a1bSJed Brown if (user->overlap) { 91c4762a1bSJed Brown DM odm = NULL; 92c4762a1bSJed Brown /* Add the level-1 overlap to refined mesh */ 93c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_OVERLAP]);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = DMPlexDistributeOverlap(*dm, 1, NULL, &odm);CHKERRQ(ierr); 95c4762a1bSJed Brown if (odm) { 96c4762a1bSJed Brown ierr = DMView(odm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 98c4762a1bSJed Brown *dm = odm; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown if (simplex){ierr = PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh");CHKERRQ(ierr);} 103c4762a1bSJed Brown else{ierr = PetscObjectSetName((PetscObject) *dm, "Tensor Product Mesh");CHKERRQ(ierr);} 104c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 106c4762a1bSJed Brown PetscFunctionReturn(0); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown int main(int argc, char **argv) 110c4762a1bSJed Brown { 111c4762a1bSJed Brown DM dm, pdm; 112c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 113c4762a1bSJed Brown PetscErrorCode ierr; 114c4762a1bSJed Brown PetscPartitioner part; 115c4762a1bSJed Brown 116c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 117c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = DMPlexDistribute(dm, user.overlap, NULL, &pdm);CHKERRQ(ierr); 122c4762a1bSJed Brown if (pdm) {ierr = DMViewFromOptions(pdm, NULL, "-pdm_view");CHKERRQ(ierr);} 123c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = DMDestroy(&pdm);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = PetscFinalize(); 126c4762a1bSJed Brown return ierr; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown /*TEST 130c4762a1bSJed Brown 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: 0 133c4762a1bSJed Brown requires: ctetgen 134c4762a1bSJed Brown args: -dim 3 -dm_refine 2 -petscpartitioner_type simple -dm_view 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: 1 137c4762a1bSJed Brown args: -dim 3 -simplex 0 -dm_refine 2 -petscpartitioner_type simple -dm_view 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: quad_0 140c4762a1bSJed Brown nsize: 2 141c4762a1bSJed Brown args: -dim 3 -simplex 0 -dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view 142c4762a1bSJed Brown 143c4762a1bSJed Brown TEST*/ 144