1c4762a1bSJed Brown static char help[] = "Orient a mesh in parallel\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown /* Domain and mesh definition */ 7c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 8c4762a1bSJed Brown PetscInt testNum; /* Labels the different test partitions */ 9c4762a1bSJed Brown } AppCtx; 10c4762a1bSJed Brown 11c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown PetscErrorCode ierr; 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscFunctionBeginUser; 16c4762a1bSJed Brown options->testPartition = PETSC_TRUE; 17c4762a1bSJed Brown options->testNum = 0; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 20c4762a1bSJed Brown ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex13.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr); 21c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-test_num", "The test partition number", "ex13.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr); 221e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 23c4762a1bSJed Brown PetscFunctionReturn(0); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 26c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 27c4762a1bSJed Brown { 28c4762a1bSJed Brown DM dmDist = NULL; 2930602db0SMatthew G. Knepley PetscBool simplex; 3030602db0SMatthew G. Knepley PetscInt dim; 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscFunctionBeginUser; 3430602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 3530602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 3630602db0SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 3730602db0SMatthew G. Knepley ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 3830602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(*dm, &simplex);CHKERRQ(ierr); 39c4762a1bSJed Brown if (user->testPartition) { 40c4762a1bSJed Brown PetscPartitioner part; 41c4762a1bSJed Brown PetscInt *sizes = NULL; 42c4762a1bSJed Brown PetscInt *points = NULL; 43c4762a1bSJed Brown PetscMPIInt rank, size; 44c4762a1bSJed Brown 45ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 46ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 47dd400576SPatrick Sanan if (rank == 0) { 4830602db0SMatthew G. Knepley if (dim == 2 && simplex && size == 2) { 49c4762a1bSJed Brown switch (user->testNum) { 50c4762a1bSJed Brown case 0: { 51c4762a1bSJed Brown PetscInt triSizes_p2[2] = {4, 4}; 52c4762a1bSJed Brown PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4}; 53c4762a1bSJed Brown 54c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscArraycpy(sizes, triSizes_p2, 2);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr); 57c4762a1bSJed Brown break;} 58c4762a1bSJed Brown case 1: { 59c4762a1bSJed Brown PetscInt triSizes_p2[2] = {6, 2}; 60c4762a1bSJed Brown PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5}; 61c4762a1bSJed Brown 62c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = PetscArraycpy(sizes, triSizes_p2, 2);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr); 65c4762a1bSJed Brown break;} 66c4762a1bSJed Brown default: 67*98921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum); 68c4762a1bSJed Brown } 6930602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 3) { 70c4762a1bSJed Brown PetscInt triSizes_p3[3] = {3, 3, 2}; 71c4762a1bSJed Brown PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5}; 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = PetscMalloc2(3, &sizes, 8, &points);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = PetscArraycpy(sizes, triSizes_p3, 3);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = PetscArraycpy(points, triPoints_p3, 8);CHKERRQ(ierr); 7630602db0SMatthew G. Knepley } else if (dim == 2 && !simplex && size == 2) { 77c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {2, 2}; 78c4762a1bSJed Brown PetscInt quadPoints_p2[4] = {2, 3, 0, 1}; 79c4762a1bSJed Brown 80c4762a1bSJed Brown ierr = PetscMalloc2(2, &sizes, 4, &points);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = PetscArraycpy(sizes, quadSizes_p2, 2);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = PetscArraycpy(points, quadPoints_p2, 4);CHKERRQ(ierr); 83c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition"); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = PetscFree2(sizes, points);CHKERRQ(ierr); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr); 91c4762a1bSJed Brown if (dmDist) { 92c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 93c4762a1bSJed Brown *dm = dmDist; 94c4762a1bSJed Brown } 9530602db0SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, simplex ? "Simplicial Mesh" : "Tensor Product Mesh");CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 97c4762a1bSJed Brown PetscFunctionReturn(0); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown static PetscErrorCode ScrambleOrientation(DM dm, AppCtx *user) 101c4762a1bSJed Brown { 102c4762a1bSJed Brown PetscInt h, cStart, cEnd, c; 103c4762a1bSJed Brown PetscErrorCode ierr; 104c4762a1bSJed Brown 105c4762a1bSJed Brown PetscFunctionBeginUser; 106c4762a1bSJed Brown ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); 108c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 109c4762a1bSJed Brown /* Could use PetscRand instead */ 110b5a892a1SMatthew G. Knepley if (c%2) {ierr = DMPlexOrientPoint(dm, c, -1);CHKERRQ(ierr);} 111c4762a1bSJed Brown } 112c4762a1bSJed Brown PetscFunctionReturn(0); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown static PetscErrorCode TestOrientation(DM dm, AppCtx *user) 116c4762a1bSJed Brown { 117c4762a1bSJed Brown PetscErrorCode ierr; 118c4762a1bSJed Brown 119c4762a1bSJed Brown PetscFunctionBeginUser; 120c4762a1bSJed Brown ierr = ScrambleOrientation(dm, user);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = DMPlexOrient(dm);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = DMViewFromOptions(dm, NULL, "-oriented_dm_view");CHKERRQ(ierr); 123c4762a1bSJed Brown PetscFunctionReturn(0); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown int main(int argc, char **argv) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown DM dm; 129c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 130c4762a1bSJed Brown PetscErrorCode ierr; 131c4762a1bSJed Brown 132c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 133c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = TestOrientation(dm, &user);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = PetscFinalize(); 138c4762a1bSJed Brown return ierr; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /*TEST 14230602db0SMatthew G. Knepley testset: 14330602db0SMatthew G. Knepley requires: triangle 14430602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -oriented_dm_view ascii::ascii_info_detail -orientation_view 14530602db0SMatthew G. Knepley 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown suffix: 0 14830602db0SMatthew G. Knepley args: -test_partition 0 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown suffix: 1 151c4762a1bSJed Brown nsize: 2 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: 2 154c4762a1bSJed Brown nsize: 2 15530602db0SMatthew G. Knepley args: -test_num 1 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: 3 158c4762a1bSJed Brown nsize: 3 15930602db0SMatthew G. Knepley args: -orientation_view_synchronized 160c4762a1bSJed Brown 161c4762a1bSJed Brown TEST*/ 162