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