xref: /petsc/src/dm/tutorials/ex14.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscFunctionBeginUser;
23*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   if (info.dim == 3) {
27c4762a1bSJed Brown     PetscScalar    ***g;
28*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,gvec,&g));
29c4762a1bSJed Brown     /* loop over ghosts */
30c4762a1bSJed Brown     for (k=info.zs; k<info.zs+info.zm; k++) {
31c4762a1bSJed Brown       for (j=info.ys; j<info.ys+info.ym; j++) {
32c4762a1bSJed Brown         for (i=info.xs; i<info.xs+info.xm; i++) {
33c4762a1bSJed Brown           g[k][j][info.dof*i+0]   = i;
34c4762a1bSJed Brown           g[k][j][info.dof*i+1]   = j;
35c4762a1bSJed Brown           g[k][j][info.dof*i+2]   = k;
36c4762a1bSJed Brown         }
37c4762a1bSJed Brown       }
38c4762a1bSJed Brown     }
39*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,gvec,&g));
40c4762a1bSJed Brown   }
41c4762a1bSJed Brown   if (info.dim == 2) {
42c4762a1bSJed Brown     PetscScalar    **g;
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,gvec,&g));
44c4762a1bSJed Brown     /* loop over ghosts */
45c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
46c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
47c4762a1bSJed Brown         for (l = 0;l<info.dof;l++) {
48c4762a1bSJed Brown           g[j][info.dof*i+0]   = i;
49c4762a1bSJed Brown           g[j][info.dof*i+1]   = j;
50c4762a1bSJed Brown           g[j][info.dof*i+2]   = rank;
51c4762a1bSJed Brown         }
52c4762a1bSJed Brown       }
53c4762a1bSJed Brown     }
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,gvec,&g));
55c4762a1bSJed Brown   }
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown int main(int argc,char **argv)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   PetscErrorCode ierr;
62c4762a1bSJed Brown   DM             da,*subda;
63c4762a1bSJed Brown   PetscInt       i,dim = 3;
643b5e53cdSSajid Ali   PetscInt       M = 25, N = 25, P = 25;
65c4762a1bSJed Brown   PetscMPIInt    size,rank;
66c4762a1bSJed Brown   Vec            v;
67c4762a1bSJed Brown   Vec            slvec,sgvec;
68c4762a1bSJed Brown   IS             *ois,*iis;
69c4762a1bSJed Brown   VecScatter     oscata;
70c4762a1bSJed Brown   VecScatter     *iscat,*oscat,*gscat;
71c4762a1bSJed Brown   DMDALocalInfo  info;
723b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_TRUE;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Create distributed array and get vectors */
78*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
79*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
80c4762a1bSJed Brown   if (dim == 2) {
81*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,&da));
82c4762a1bSJed Brown   } else if (dim == 3) {
83*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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));
84c4762a1bSJed Brown   }
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
88c4762a1bSJed Brown 
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDomainDecomposition(da,NULL,NULL,&iis,&ois,&subda));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDomainDecompositionScatters(da,1,subda,&iscat,&oscat,&gscat));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   {
93c4762a1bSJed Brown     DMDALocalInfo subinfo;
94c4762a1bSJed Brown     MatStencil    lower,upper;
953b5e53cdSSajid Ali     IS            patchis;
96c4762a1bSJed Brown     Vec           smallvec;
97c4762a1bSJed Brown     Vec           largevec;
98c4762a1bSJed Brown     VecScatter    patchscat;
99c4762a1bSJed Brown 
100*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetLocalInfo(subda[0],&subinfo));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown     lower.i = info.xs;
103c4762a1bSJed Brown     lower.j = info.ys;
104c4762a1bSJed Brown     lower.k = info.zs;
105c4762a1bSJed Brown     upper.i = info.xs+info.xm;
106c4762a1bSJed Brown     upper.j = info.ys+info.ym;
107c4762a1bSJed Brown     upper.k = info.zs+info.zm;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown     /* test the patch IS as a thing to scatter to/from */
110*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreatePatchIS(da,&lower,&upper,&patchis,patchis_offproc));
111*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(da,&largevec));
112c4762a1bSJed Brown 
113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreate(PETSC_COMM_SELF,&smallvec));
114*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetSizes(smallvec,info.dof*(upper.i - lower.i)*(upper.j - lower.j)*(upper.k - lower.k),PETSC_DECIDE));
115*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetFromOptions(smallvec));
116*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(smallvec,NULL,largevec,patchis,&patchscat));
117c4762a1bSJed Brown 
118*5f80ce2aSJacob Faibussowitsch     CHKERRQ(FillLocalSubdomain(subda[0],smallvec));
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(largevec,0));
120c4762a1bSJed Brown 
121*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD));
122*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD));
123*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(patchis,PETSC_VIEWER_STDOUT_WORLD));
124*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterView(patchscat,PETSC_VIEWER_STDOUT_WORLD));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown     for (i = 0; i < size; i++) {
127c4762a1bSJed Brown       if (i == rank) {
128*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecView(smallvec,PETSC_VIEWER_STDOUT_SELF));
129c4762a1bSJed Brown       }
130*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
131c4762a1bSJed Brown     }
132c4762a1bSJed Brown 
133*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
134*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(largevec,PETSC_VIEWER_STDOUT_WORLD));
135c4762a1bSJed Brown 
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&smallvec));
137*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreGlobalVector(da,&largevec));
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&patchis));
139*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&patchscat));
140c4762a1bSJed Brown   }
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* view the various parts */
143c4762a1bSJed Brown   {
144c4762a1bSJed Brown     for (i = 0; i < size; i++) {
145c4762a1bSJed Brown       if (i == rank) {
146*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i));
147*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMView(subda[0],PETSC_VIEWER_STDOUT_SELF));
148c4762a1bSJed Brown       }
149*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
150c4762a1bSJed Brown     }
151c4762a1bSJed Brown 
152*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(subda[0],&slvec));
153*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(subda[0],&sgvec));
154*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(da,&v));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown     /* test filling outer between the big DM and the small ones with the IS scatter*/
157*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(v,ois[0],sgvec,NULL,&oscata));
158c4762a1bSJed Brown 
159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(FillLocalSubdomain(subda[0],sgvec));
160c4762a1bSJed Brown 
161*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE));
162*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown     /* test the local-to-local scatter */
165c4762a1bSJed Brown 
166c4762a1bSJed Brown     /* fill up the local subdomain and then add them together */
167*5f80ce2aSJacob Faibussowitsch     CHKERRQ(FillLocalSubdomain(da,v));
168c4762a1bSJed Brown 
169*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD));
170*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD));
171c4762a1bSJed Brown 
172*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     /* test ghost scattering backwards */
175c4762a1bSJed Brown 
176*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(v,0));
177c4762a1bSJed Brown 
178*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE));
179*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE));
180c4762a1bSJed Brown 
181*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown     /* test overlap scattering backwards */
184c4762a1bSJed Brown 
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalBegin(subda[0],slvec,ADD_VALUES,sgvec));
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalEnd(subda[0],slvec,ADD_VALUES,sgvec));
187c4762a1bSJed Brown 
188*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(v,0));
189c4762a1bSJed Brown 
190*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
191*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
192c4762a1bSJed Brown 
193*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown     /* test interior scattering backwards */
196c4762a1bSJed Brown 
197*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(v,0));
198c4762a1bSJed Brown 
199*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
200*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE));
201c4762a1bSJed Brown 
202*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown     /* test matrix allocation */
205c4762a1bSJed Brown     for (i = 0; i < size; i++) {
206c4762a1bSJed Brown       if (i == rank) {
207c4762a1bSJed Brown         Mat m;
208*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i));
209*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSetMatType(subda[0],MATAIJ));
210*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMCreateMatrix(subda[0],&m));
211*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatView(m,PETSC_VIEWER_STDOUT_SELF));
212*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatDestroy(&m));
213c4762a1bSJed Brown       }
214*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
215c4762a1bSJed Brown     }
216*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(subda[0],&slvec));
217*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreGlobalVector(subda[0],&sgvec));
218*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreGlobalVector(da,&v));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&subda[0]));
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ois[0]));
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iis[0]));
224c4762a1bSJed Brown 
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&iscat[0]));
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&oscat[0]));
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&gscat[0]));
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&oscata));
229c4762a1bSJed Brown 
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(iscat));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(oscat));
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(gscat));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(oscata));
234c4762a1bSJed Brown 
235*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(subda));
236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ois));
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(iis));
238c4762a1bSJed Brown 
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
240c4762a1bSJed Brown   ierr = PetscFinalize();
241c4762a1bSJed Brown   return ierr;
242c4762a1bSJed Brown }
243