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