xref: /petsc/src/dm/tests/ex16.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 static char help[] = "Tests DMComposite routines.\n\n";
3 
4 #include <petscdmredundant.h>
5 #include <petscdm.h>
6 #include <petscdmda.h>
7 #include <petscdmcomposite.h>
8 #include <petscpf.h>
9 
10 int main(int argc,char **argv)
11 {
12   PetscInt               nredundant1 = 5,nredundant2 = 2,i;
13   ISLocalToGlobalMapping *ltog;
14   PetscMPIInt            rank,size;
15   DM                     packer;
16   Vec                    global,local1,local2,redundant1,redundant2;
17   PF                     pf;
18   DM                     da1,da2,dmred1,dmred2;
19   PetscScalar            *redundant1a,*redundant2a;
20   PetscViewer            sviewer;
21   PetscBool              gather_add = PETSC_FALSE;
22 
23   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
24   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
25   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
26 
27   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-gather_add",&gather_add,NULL));
28 
29   CHKERRQ(DMCompositeCreate(PETSC_COMM_WORLD,&packer));
30 
31   CHKERRQ(DMRedundantCreate(PETSC_COMM_WORLD,0,nredundant1,&dmred1));
32   CHKERRQ(DMCreateLocalVector(dmred1,&redundant1));
33   CHKERRQ(DMCompositeAddDM(packer,dmred1));
34 
35   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da1));
36   CHKERRQ(DMSetFromOptions(da1));
37   CHKERRQ(DMSetUp(da1));
38   CHKERRQ(DMCreateLocalVector(da1,&local1));
39   CHKERRQ(DMCompositeAddDM(packer,da1));
40 
41   CHKERRQ(DMRedundantCreate(PETSC_COMM_WORLD,1%size,nredundant2,&dmred2));
42   CHKERRQ(DMCreateLocalVector(dmred2,&redundant2));
43   CHKERRQ(DMCompositeAddDM(packer,dmred2));
44 
45   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,6,1,1,NULL,&da2));
46   CHKERRQ(DMSetFromOptions(da2));
47   CHKERRQ(DMSetUp(da2));
48   CHKERRQ(DMCreateLocalVector(da2,&local2));
49   CHKERRQ(DMCompositeAddDM(packer,da2));
50 
51   CHKERRQ(DMCreateGlobalVector(packer,&global));
52   CHKERRQ(PFCreate(PETSC_COMM_WORLD,1,1,&pf));
53   CHKERRQ(PFSetType(pf,PFIDENTITY,NULL));
54   CHKERRQ(PFApplyVec(pf,NULL,global));
55   CHKERRQ(PFDestroy(&pf));
56   CHKERRQ(VecView(global,PETSC_VIEWER_STDOUT_WORLD));
57 
58   CHKERRQ(DMCompositeScatter(packer,global,redundant1,local1,redundant2,local2));
59   CHKERRQ(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
60   CHKERRQ(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of redundant1 vector\n",rank));
61   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
62   CHKERRQ(VecView(redundant1,sviewer));
63   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
64   CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
65   CHKERRQ(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of da1 vector\n",rank));
66   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
67   CHKERRQ(VecView(local1,sviewer));
68   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
69   CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
70   CHKERRQ(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of redundant2 vector\n",rank));
71   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
72   CHKERRQ(VecView(redundant2,sviewer));
73   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
74   CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
75   CHKERRQ(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of da2 vector\n",rank));
76   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
77   CHKERRQ(VecView(local2,sviewer));
78   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
79   CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
80   CHKERRQ(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
81 
82   CHKERRQ(VecGetArray(redundant1,&redundant1a));
83   CHKERRQ(VecGetArray(redundant2,&redundant2a));
84   for (i=0; i<nredundant1; i++) redundant1a[i] = (rank+2)*i;
85   for (i=0; i<nredundant2; i++) redundant2a[i] = (rank+10)*i;
86   CHKERRQ(VecRestoreArray(redundant1,&redundant1a));
87   CHKERRQ(VecRestoreArray(redundant2,&redundant2a));
88 
89   CHKERRQ(DMCompositeGather(packer,gather_add ? ADD_VALUES : INSERT_VALUES,global,redundant1,local1,redundant2,local2));
90   CHKERRQ(VecView(global,PETSC_VIEWER_STDOUT_WORLD));
91 
92   /* get the global numbering for each subvector element */
93   CHKERRQ(DMCompositeGetISLocalToGlobalMappings(packer,&ltog));
94 
95   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of redundant1 vector\n"));
96   CHKERRQ(ISLocalToGlobalMappingView(ltog[0],PETSC_VIEWER_STDOUT_WORLD));
97   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of local1 vector\n"));
98   CHKERRQ(ISLocalToGlobalMappingView(ltog[1],PETSC_VIEWER_STDOUT_WORLD));
99   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of redundant2 vector\n"));
100   CHKERRQ(ISLocalToGlobalMappingView(ltog[2],PETSC_VIEWER_STDOUT_WORLD));
101   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of local2 vector\n"));
102   CHKERRQ(ISLocalToGlobalMappingView(ltog[3],PETSC_VIEWER_STDOUT_WORLD));
103 
104   for (i=0; i<4; i++) CHKERRQ(ISLocalToGlobalMappingDestroy(&ltog[i]));
105   CHKERRQ(PetscFree(ltog));
106 
107   CHKERRQ(DMDestroy(&da1));
108   CHKERRQ(DMDestroy(&dmred1));
109   CHKERRQ(DMDestroy(&dmred2));
110   CHKERRQ(DMDestroy(&da2));
111   CHKERRQ(VecDestroy(&redundant1));
112   CHKERRQ(VecDestroy(&redundant2));
113   CHKERRQ(VecDestroy(&local1));
114   CHKERRQ(VecDestroy(&local2));
115   CHKERRQ(VecDestroy(&global));
116   CHKERRQ(DMDestroy(&packer));
117   CHKERRQ(PetscFinalize());
118   return 0;
119 }
120 
121 /*TEST
122 
123    build:
124       requires: !complex
125 
126    test:
127       nsize: 3
128 
129    test:
130       suffix: 2
131       nsize: 3
132       args: -gather_add
133 
134 TEST*/
135