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 11*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 12c4762a1bSJed Brown PetscFunctionBeginUser; 13c4762a1bSJed Brown options->testPartition = PETSC_TRUE; 14c4762a1bSJed Brown options->testNum = 0; 15c4762a1bSJed Brown 16d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex13.c", options->testPartition, &options->testPartition, NULL)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-test_num", "The test partition number", "ex13.c", options->testNum, &options->testNum, NULL, 0)); 19d0609cedSBarry Smith PetscOptionsEnd(); 20c4762a1bSJed Brown PetscFunctionReturn(0); 21c4762a1bSJed Brown } 22c4762a1bSJed Brown 23*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 24c4762a1bSJed Brown DM dmDist = NULL; 2530602db0SMatthew G. Knepley PetscBool simplex; 2630602db0SMatthew G. Knepley PetscInt dim; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscFunctionBeginUser; 299566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 309566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 319566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 329566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 339566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 349566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(*dm, &simplex)); 35c4762a1bSJed Brown if (user->testPartition) { 36c4762a1bSJed Brown PetscPartitioner part; 37c4762a1bSJed Brown PetscInt *sizes = NULL; 38c4762a1bSJed Brown PetscInt *points = NULL; 39c4762a1bSJed Brown PetscMPIInt rank, size; 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 43dd400576SPatrick Sanan if (rank == 0) { 4430602db0SMatthew G. Knepley if (dim == 2 && simplex && size == 2) { 45c4762a1bSJed Brown switch (user->testNum) { 46c4762a1bSJed Brown case 0: { 47c4762a1bSJed Brown PetscInt triSizes_p2[2] = {4, 4}; 48c4762a1bSJed Brown PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4}; 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 8, &points)); 519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, triSizes_p2, 2)); 529566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(points, triPoints_p2, 8)); 53*9371c9d4SSatish Balay break; 54*9371c9d4SSatish Balay } 55c4762a1bSJed Brown case 1: { 56c4762a1bSJed Brown PetscInt triSizes_p2[2] = {6, 2}; 57c4762a1bSJed Brown PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5}; 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 8, &points)); 609566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, triSizes_p2, 2)); 619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(points, triPoints_p2, 8)); 62*9371c9d4SSatish Balay break; 63*9371c9d4SSatish Balay } 64*9371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum); 65c4762a1bSJed Brown } 6630602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 3) { 67c4762a1bSJed Brown PetscInt triSizes_p3[3] = {3, 3, 2}; 68c4762a1bSJed Brown PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5}; 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(3, &sizes, 8, &points)); 719566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, triSizes_p3, 3)); 729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(points, triPoints_p3, 8)); 7330602db0SMatthew G. Knepley } else if (dim == 2 && !simplex && size == 2) { 74c4762a1bSJed Brown PetscInt quadSizes_p2[2] = {2, 2}; 75c4762a1bSJed Brown PetscInt quadPoints_p2[4] = {2, 3, 0, 1}; 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2, &sizes, 4, &points)); 789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2)); 799566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(points, quadPoints_p2, 4)); 80c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition"); 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 839566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 849566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points)); 859566063dSJacob Faibussowitsch PetscCall(PetscFree2(sizes, points)); 86c4762a1bSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, 0, NULL, &dmDist)); 88c4762a1bSJed Brown if (dmDist) { 899566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 90c4762a1bSJed Brown *dm = dmDist; 91c4762a1bSJed Brown } 929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, simplex ? "Simplicial Mesh" : "Tensor Product Mesh")); 939566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 94c4762a1bSJed Brown PetscFunctionReturn(0); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97*9371c9d4SSatish Balay static PetscErrorCode ScrambleOrientation(DM dm, AppCtx *user) { 98c4762a1bSJed Brown PetscInt h, cStart, cEnd, c; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBeginUser; 1019566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &h)); 1029566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd)); 103c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 104c4762a1bSJed Brown /* Could use PetscRand instead */ 1059566063dSJacob Faibussowitsch if (c % 2) PetscCall(DMPlexOrientPoint(dm, c, -1)); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown PetscFunctionReturn(0); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110*9371c9d4SSatish Balay static PetscErrorCode TestOrientation(DM dm, AppCtx *user) { 111c4762a1bSJed Brown PetscFunctionBeginUser; 1129566063dSJacob Faibussowitsch PetscCall(ScrambleOrientation(dm, user)); 1139566063dSJacob Faibussowitsch PetscCall(DMPlexOrient(dm)); 1149566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-oriented_dm_view")); 115c4762a1bSJed Brown PetscFunctionReturn(0); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118*9371c9d4SSatish Balay int main(int argc, char **argv) { 119c4762a1bSJed Brown DM dm; 120c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 121c4762a1bSJed Brown 122327415f7SBarry Smith PetscFunctionBeginUser; 1239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1249566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1259566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1269566063dSJacob Faibussowitsch PetscCall(TestOrientation(dm, &user)); 1279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1289566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 129b122ec5aSJacob Faibussowitsch return 0; 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /*TEST 13330602db0SMatthew G. Knepley testset: 13430602db0SMatthew G. Knepley requires: triangle 13530602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -oriented_dm_view ascii::ascii_info_detail -orientation_view 13630602db0SMatthew G. Knepley 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: 0 13930602db0SMatthew G. Knepley args: -test_partition 0 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: 1 142c4762a1bSJed Brown nsize: 2 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: 2 145c4762a1bSJed Brown nsize: 2 14630602db0SMatthew G. Knepley args: -test_num 1 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: 3 149c4762a1bSJed Brown nsize: 3 15030602db0SMatthew G. Knepley args: -orientation_view_synchronized 151c4762a1bSJed Brown 152c4762a1bSJed Brown TEST*/ 153