xref: /petsc/src/dm/tests/ex51.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 
2 static char help[] = "Tests DMDAGlobalToNaturalAllCreate() using contour plotting for 2d DMDAs.\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 #include <petscdraw.h>
7 
8 int main(int argc, char **argv)
9 {
10   PetscInt        i, j, M = 10, N = 8, m = PETSC_DECIDE, n = PETSC_DECIDE;
11   PetscMPIInt     rank;
12   PetscBool       flg = PETSC_FALSE;
13   DM              da;
14   PetscViewer     viewer;
15   Vec             localall, global;
16   PetscScalar     value, *vlocal;
17   DMBoundaryType  bx = DM_BOUNDARY_NONE, by = DM_BOUNDARY_NONE;
18   DMDAStencilType stype = DMDA_STENCIL_BOX;
19   VecScatter      tolocalall, fromlocalall;
20   PetscInt        start, end;
21 
22   PetscFunctionBeginUser;
23   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
24   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "", 300, 0, 300, 300, &viewer));
25 
26   /* Read options */
27   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
28   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
29   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
30   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
31   PetscCall(PetscOptionsGetBool(NULL, NULL, "-star_stencil", &flg, NULL));
32   if (flg) stype = DMDA_STENCIL_STAR;
33 
34   /* Create distributed array and get vectors */
35   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, m, n, 1, 1, NULL, NULL, &da));
36   PetscCall(DMSetFromOptions(da));
37   PetscCall(DMSetUp(da));
38 
39   PetscCall(DMCreateGlobalVector(da, &global));
40   PetscCall(VecCreateSeq(PETSC_COMM_SELF, M * N, &localall));
41 
42   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
43   PetscCall(VecGetOwnershipRange(global, &start, &end));
44   for (i = start; i < end; i++) {
45     value = 5.0 * rank;
46     PetscCall(VecSetValues(global, 1, &i, &value, INSERT_VALUES));
47   }
48   PetscCall(VecView(global, viewer));
49 
50   /*
51      Create Scatter from global DMDA parallel vector to local vector that
52    contains all entries
53   */
54   PetscCall(DMDAGlobalToNaturalAllCreate(da, &tolocalall));
55   PetscCall(DMDANaturalAllToGlobalCreate(da, &fromlocalall));
56 
57   PetscCall(VecScatterBegin(tolocalall, global, localall, INSERT_VALUES, SCATTER_FORWARD));
58   PetscCall(VecScatterEnd(tolocalall, global, localall, INSERT_VALUES, SCATTER_FORWARD));
59 
60   PetscCall(VecGetArray(localall, &vlocal));
61   for (j = 0; j < N; j++) {
62     for (i = 0; i < M; i++) *vlocal++ += i + j * M;
63   }
64   PetscCall(VecRestoreArray(localall, &vlocal));
65 
66   /* scatter back to global vector */
67   PetscCall(VecScatterBegin(fromlocalall, localall, global, INSERT_VALUES, SCATTER_FORWARD));
68   PetscCall(VecScatterEnd(fromlocalall, localall, global, INSERT_VALUES, SCATTER_FORWARD));
69 
70   PetscCall(VecView(global, viewer));
71 
72   /* Free memory */
73   PetscCall(VecScatterDestroy(&tolocalall));
74   PetscCall(VecScatterDestroy(&fromlocalall));
75   PetscCall(PetscViewerDestroy(&viewer));
76   PetscCall(VecDestroy(&localall));
77   PetscCall(VecDestroy(&global));
78   PetscCall(DMDestroy(&da));
79   PetscCall(PetscFinalize());
80   return 0;
81 }
82 
83 /*TEST
84 
85    build:
86      requires: !complex
87 
88    test:
89       nsize: 3
90 
91 TEST*/
92