xref: /petsc/src/dm/tests/ex44.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1 static char help[] = "Tests various DMComposite routines.\n\n";
2 
3 #include <petscdm.h>
4 #include <petscdmda.h>
5 #include <petscdmcomposite.h>
6 
main(int argc,char ** argv)7 int main(int argc, char **argv)
8 {
9   PetscMPIInt rank;
10   DM          da1, da2, packer;
11   Vec         local, global, globals[2], buffer;
12   PetscScalar value, *array;
13   PetscInt    n;
14   PetscViewer viewer;
15 
16   PetscFunctionBeginUser;
17   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
18 
19   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &packer));
20   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 8, 1, 1, NULL, &da1));
21   PetscCall(DMSetFromOptions(da1));
22   PetscCall(DMSetUp(da1));
23   PetscCall(DMCompositeAddDM(packer, da1));
24   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 6, 1, 1, NULL, &da2));
25   PetscCall(DMSetFromOptions(da2));
26   PetscCall(DMSetUp(da2));
27   PetscCall(DMCompositeAddDM(packer, da2));
28 
29   PetscCall(DMCreateGlobalVector(packer, &global));
30   PetscCall(DMCreateLocalVector(packer, &local));
31   PetscCall(DMCreateLocalVector(packer, &buffer));
32 
33   PetscCall(DMCompositeGetAccessArray(packer, global, 2, NULL, globals));
34   value = 1;
35   PetscCall(VecSet(globals[0], value));
36   value = -1;
37   PetscCall(VecSet(globals[1], value));
38   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
39   value = rank + 1;
40   PetscCall(VecGetLocalSize(globals[0], &n));
41   PetscCall(VecGetArray(globals[0], &array));
42   for (PetscInt i = 0; i < n; i++) array[i] *= value;
43   PetscCall(VecRestoreArray(globals[0], &array));
44   PetscCall(VecGetLocalSize(globals[1], &n));
45   PetscCall(VecGetArray(globals[1], &array));
46   for (PetscInt i = 0; i < n; i++) array[i] *= value;
47   PetscCall(VecRestoreArray(globals[1], &array));
48   PetscCall(DMCompositeRestoreAccessArray(packer, global, 2, NULL, globals));
49 
50   /* Test GlobalToLocal in insert mode */
51   PetscCall(DMGlobalToLocalBegin(packer, global, INSERT_VALUES, local));
52   PetscCall(DMGlobalToLocalEnd(packer, global, INSERT_VALUES, local));
53 
54   PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
55   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
56   PetscCall(PetscViewerASCIIPrintf(viewer, "\nLocal Vector: processor %d\n", rank));
57   PetscCall(VecView(local, viewer));
58   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
59   PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
60 
61   /* Test LocalToGlobal in insert mode */
62   PetscCall(DMLocalToGlobalBegin(packer, local, INSERT_VALUES, global));
63   PetscCall(DMLocalToGlobalEnd(packer, local, INSERT_VALUES, global));
64 
65   PetscCall(VecView(global, PETSC_VIEWER_STDOUT_WORLD));
66 
67   /* Test LocalToLocal in insert mode */
68   PetscCall(DMLocalToLocalBegin(packer, local, INSERT_VALUES, buffer));
69   PetscCall(DMLocalToLocalEnd(packer, local, INSERT_VALUES, buffer));
70 
71   PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
72   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
73   PetscCall(PetscViewerASCIIPrintf(viewer, "\nLocal Vector: processor %d\n", rank));
74   PetscCall(VecView(buffer, viewer));
75   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &viewer));
76   PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
77 
78   PetscCall(VecDestroy(&buffer));
79   PetscCall(VecDestroy(&local));
80   PetscCall(VecDestroy(&global));
81   PetscCall(DMDestroy(&packer));
82   PetscCall(DMDestroy(&da2));
83   PetscCall(DMDestroy(&da1));
84 
85   PetscCall(PetscFinalize());
86   return 0;
87 }
88 
89 /*TEST
90 
91    test:
92       nsize: 3
93 
94 TEST*/
95