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