xref: /petsc/src/dm/tutorials/ex14.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests DMCreateDomainDecomposition.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown Use the options
6c4762a1bSJed Brown      -da_grid_x <nx> - number of grid points in x direction, if M < 0
7c4762a1bSJed Brown      -da_grid_y <ny> - number of grid points in y direction, if N < 0
8c4762a1bSJed Brown      -da_processors_x <MX> number of processors in x directio
9c4762a1bSJed Brown      -da_processors_y <MY> number of processors in x direction
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petscdm.h>
13c4762a1bSJed Brown #include <petscdmda.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown PetscErrorCode FillLocalSubdomain(DM da, Vec gvec)
16c4762a1bSJed Brown {
17c4762a1bSJed Brown   DMDALocalInfo  info;
18c4762a1bSJed Brown   PetscMPIInt    rank;
19c4762a1bSJed Brown   PetscInt       i,j,k,l;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
239566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   if (info.dim == 3) {
26c4762a1bSJed Brown     PetscScalar    ***g;
279566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da,gvec,&g));
28c4762a1bSJed Brown     /* loop over ghosts */
29c4762a1bSJed Brown     for (k=info.zs; k<info.zs+info.zm; k++) {
30c4762a1bSJed Brown       for (j=info.ys; j<info.ys+info.ym; j++) {
31c4762a1bSJed Brown         for (i=info.xs; i<info.xs+info.xm; i++) {
32c4762a1bSJed Brown           g[k][j][info.dof*i+0]   = i;
33c4762a1bSJed Brown           g[k][j][info.dof*i+1]   = j;
34c4762a1bSJed Brown           g[k][j][info.dof*i+2]   = k;
35c4762a1bSJed Brown         }
36c4762a1bSJed Brown       }
37c4762a1bSJed Brown     }
389566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da,gvec,&g));
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown   if (info.dim == 2) {
41c4762a1bSJed Brown     PetscScalar    **g;
429566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da,gvec,&g));
43c4762a1bSJed Brown     /* loop over ghosts */
44c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
45c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
46c4762a1bSJed Brown         for (l = 0;l<info.dof;l++) {
47c4762a1bSJed Brown           g[j][info.dof*i+0]   = i;
48c4762a1bSJed Brown           g[j][info.dof*i+1]   = j;
49c4762a1bSJed Brown           g[j][info.dof*i+2]   = rank;
50c4762a1bSJed Brown         }
51c4762a1bSJed Brown       }
52c4762a1bSJed Brown     }
539566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da,gvec,&g));
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown   PetscFunctionReturn(0);
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown int main(int argc,char **argv)
59c4762a1bSJed Brown {
60c4762a1bSJed Brown   DM             da,*subda;
61c4762a1bSJed Brown   PetscInt       i,dim = 3;
623b5e53cdSSajid Ali   PetscInt       M = 25, N = 25, P = 25;
63c4762a1bSJed Brown   PetscMPIInt    size,rank;
64c4762a1bSJed Brown   Vec            v;
65c4762a1bSJed Brown   Vec            slvec,sgvec;
66c4762a1bSJed Brown   IS             *ois,*iis;
67c4762a1bSJed Brown   VecScatter     oscata;
68c4762a1bSJed Brown   VecScatter     *iscat,*oscat,*gscat;
69c4762a1bSJed Brown   DMDALocalInfo  info;
703b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_TRUE;
71c4762a1bSJed Brown 
72*327415f7SBarry Smith   PetscFunctionBeginUser;
739566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Create distributed array and get vectors */
779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
79c4762a1bSJed Brown   if (dim == 2) {
809566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,&da));
81c4762a1bSJed Brown   } else if (dim == 3) {
829566063dSJacob Faibussowitsch     PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,NULL,&da));
83c4762a1bSJed Brown   }
849566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
859566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
869566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(DMCreateDomainDecomposition(da,NULL,NULL,&iis,&ois,&subda));
899566063dSJacob Faibussowitsch   PetscCall(DMCreateDomainDecompositionScatters(da,1,subda,&iscat,&oscat,&gscat));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   {
92c4762a1bSJed Brown     DMDALocalInfo subinfo;
93c4762a1bSJed Brown     MatStencil    lower,upper;
943b5e53cdSSajid Ali     IS            patchis;
95c4762a1bSJed Brown     Vec           smallvec;
96c4762a1bSJed Brown     Vec           largevec;
97c4762a1bSJed Brown     VecScatter    patchscat;
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(subda[0],&subinfo));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown     lower.i = info.xs;
102c4762a1bSJed Brown     lower.j = info.ys;
103c4762a1bSJed Brown     lower.k = info.zs;
104c4762a1bSJed Brown     upper.i = info.xs+info.xm;
105c4762a1bSJed Brown     upper.j = info.ys+info.ym;
106c4762a1bSJed Brown     upper.k = info.zs+info.zm;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown     /* test the patch IS as a thing to scatter to/from */
1099566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(da,&lower,&upper,&patchis,patchis_offproc));
1109566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(da,&largevec));
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF,&smallvec));
1139566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(smallvec,info.dof*(upper.i - lower.i)*(upper.j - lower.j)*(upper.k - lower.k),PETSC_DECIDE));
1149566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(smallvec));
1159566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(smallvec,NULL,largevec,patchis,&patchscat));
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch     PetscCall(FillLocalSubdomain(subda[0],smallvec));
1189566063dSJacob Faibussowitsch     PetscCall(VecSet(largevec,0));
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD));
1219566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD));
1229566063dSJacob Faibussowitsch     PetscCall(ISView(patchis,PETSC_VIEWER_STDOUT_WORLD));
1239566063dSJacob Faibussowitsch     PetscCall(VecScatterView(patchscat,PETSC_VIEWER_STDOUT_WORLD));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown     for (i = 0; i < size; i++) {
126c4762a1bSJed Brown       if (i == rank) {
1279566063dSJacob Faibussowitsch         PetscCall(VecView(smallvec,PETSC_VIEWER_STDOUT_SELF));
128c4762a1bSJed Brown       }
1299566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
130c4762a1bSJed Brown     }
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
1339566063dSJacob Faibussowitsch     PetscCall(VecView(largevec,PETSC_VIEWER_STDOUT_WORLD));
134c4762a1bSJed Brown 
1359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&smallvec));
1369566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(da,&largevec));
1379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&patchis));
1389566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&patchscat));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* view the various parts */
142c4762a1bSJed Brown   {
143c4762a1bSJed Brown     for (i = 0; i < size; i++) {
144c4762a1bSJed Brown       if (i == rank) {
1459566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i));
1469566063dSJacob Faibussowitsch         PetscCall(DMView(subda[0],PETSC_VIEWER_STDOUT_SELF));
147c4762a1bSJed Brown       }
1489566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
149c4762a1bSJed Brown     }
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(subda[0],&slvec));
1529566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(subda[0],&sgvec));
1539566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(da,&v));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown     /* test filling outer between the big DM and the small ones with the IS scatter*/
1569566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(v,ois[0],sgvec,NULL,&oscata));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch     PetscCall(FillLocalSubdomain(subda[0],sgvec));
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE));
1619566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown     /* test the local-to-local scatter */
164c4762a1bSJed Brown 
165c4762a1bSJed Brown     /* fill up the local subdomain and then add them together */
1669566063dSJacob Faibussowitsch     PetscCall(FillLocalSubdomain(da,v));
167c4762a1bSJed Brown 
1689566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD));
1699566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch     PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     /* test ghost scattering backwards */
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch     PetscCall(VecSet(v,0));
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE));
1789566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE));
179c4762a1bSJed Brown 
1809566063dSJacob Faibussowitsch     PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown     /* test overlap scattering backwards */
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(subda[0],slvec,ADD_VALUES,sgvec));
1859566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(subda[0],slvec,ADD_VALUES,sgvec));
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch     PetscCall(VecSet(v,0));
188c4762a1bSJed Brown 
1899566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
1909566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
191c4762a1bSJed Brown 
1929566063dSJacob Faibussowitsch     PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown     /* test interior scattering backwards */
195c4762a1bSJed Brown 
1969566063dSJacob Faibussowitsch     PetscCall(VecSet(v,0));
197c4762a1bSJed Brown 
1989566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
1999566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
200c4762a1bSJed Brown 
2019566063dSJacob Faibussowitsch     PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown     /* test matrix allocation */
204c4762a1bSJed Brown     for (i = 0; i < size; i++) {
205c4762a1bSJed Brown       if (i == rank) {
206c4762a1bSJed Brown         Mat m;
2079566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i));
2089566063dSJacob Faibussowitsch         PetscCall(DMSetMatType(subda[0],MATAIJ));
2099566063dSJacob Faibussowitsch         PetscCall(DMCreateMatrix(subda[0],&m));
2109566063dSJacob Faibussowitsch         PetscCall(MatView(m,PETSC_VIEWER_STDOUT_SELF));
2119566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&m));
212c4762a1bSJed Brown       }
2139566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
214c4762a1bSJed Brown     }
2159566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(subda[0],&slvec));
2169566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(subda[0],&sgvec));
2179566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(da,&v));
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown 
2209566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&subda[0]));
2219566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ois[0]));
2229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iis[0]));
223c4762a1bSJed Brown 
2249566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&iscat[0]));
2259566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&oscat[0]));
2269566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&gscat[0]));
2279566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&oscata));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(iscat));
2309566063dSJacob Faibussowitsch   PetscCall(PetscFree(oscat));
2319566063dSJacob Faibussowitsch   PetscCall(PetscFree(gscat));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFree(oscata));
233c4762a1bSJed Brown 
2349566063dSJacob Faibussowitsch   PetscCall(PetscFree(subda));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFree(ois));
2369566063dSJacob Faibussowitsch   PetscCall(PetscFree(iis));
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2399566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
240b122ec5aSJacob Faibussowitsch   return 0;
241c4762a1bSJed Brown }
242