1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests DMDAGlobalToNaturalAllCreate() using contour plotting for 2d DMDAs.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdraw.h> 7c4762a1bSJed Brown 89371c9d4SSatish Balay int main(int argc, char **argv) { 9c4762a1bSJed Brown PetscInt i, j, M = 10, N = 8, m = PETSC_DECIDE, n = PETSC_DECIDE; 10c4762a1bSJed Brown PetscMPIInt rank; 11c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 12c4762a1bSJed Brown DM da; 13c4762a1bSJed Brown PetscViewer viewer; 14c4762a1bSJed Brown Vec localall, global; 15c4762a1bSJed Brown PetscScalar value, *vlocal; 16c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_NONE, by = DM_BOUNDARY_NONE; 17c4762a1bSJed Brown DMDAStencilType stype = DMDA_STENCIL_BOX; 18c4762a1bSJed Brown VecScatter tolocalall, fromlocalall; 19c4762a1bSJed Brown PetscInt start, end; 20c4762a1bSJed Brown 21327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 239566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "", 300, 0, 300, 300, &viewer)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Read options */ 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-star_stencil", &flg, NULL)); 31c4762a1bSJed Brown if (flg) stype = DMDA_STENCIL_STAR; 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Create distributed array and get vectors */ 349566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, m, n, 1, 1, NULL, NULL, &da)); 359566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 369566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &global)); 399566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, M * N, &localall)); 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 429566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(global, &start, &end)); 43c4762a1bSJed Brown for (i = start; i < end; i++) { 44c4762a1bSJed Brown value = 5.0 * rank; 459566063dSJacob Faibussowitsch PetscCall(VecSetValues(global, 1, &i, &value, INSERT_VALUES)); 46c4762a1bSJed Brown } 479566063dSJacob Faibussowitsch PetscCall(VecView(global, viewer)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Create Scatter from global DMDA parallel vector to local vector that 51c4762a1bSJed Brown contains all entries 52c4762a1bSJed Brown */ 539566063dSJacob Faibussowitsch PetscCall(DMDAGlobalToNaturalAllCreate(da, &tolocalall)); 549566063dSJacob Faibussowitsch PetscCall(DMDANaturalAllToGlobalCreate(da, &fromlocalall)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(tolocalall, global, localall, INSERT_VALUES, SCATTER_FORWARD)); 579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(tolocalall, global, localall, INSERT_VALUES, SCATTER_FORWARD)); 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(VecGetArray(localall, &vlocal)); 60c4762a1bSJed Brown for (j = 0; j < N; j++) { 61*ad540459SPierre Jolivet for (i = 0; i < M; i++) *vlocal++ += i + j * M; 62c4762a1bSJed Brown } 639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localall, &vlocal)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* scatter back to global vector */ 669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(fromlocalall, localall, global, INSERT_VALUES, SCATTER_FORWARD)); 679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(fromlocalall, localall, global, INSERT_VALUES, SCATTER_FORWARD)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(VecView(global, viewer)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Free memory */ 729566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&tolocalall)); 739566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&fromlocalall)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&localall)); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 779566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 789566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 79b122ec5aSJacob Faibussowitsch return 0; 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /*TEST 83c4762a1bSJed Brown 84c4762a1bSJed Brown build: 85c4762a1bSJed Brown requires: !complex 86c4762a1bSJed Brown 87c4762a1bSJed Brown test: 88c4762a1bSJed Brown nsize: 3 89c4762a1bSJed Brown 90c4762a1bSJed Brown TEST*/ 91