xref: /petsc/src/dm/tutorials/ex14.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Tests DMCreateDomainDecomposition.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown Use the options
5c4762a1bSJed Brown      -da_grid_x <nx> - number of grid points in x direction, if M < 0
6c4762a1bSJed Brown      -da_grid_y <ny> - number of grid points in y direction, if N < 0
7c4762a1bSJed Brown      -da_processors_x <MX> number of processors in x directio
8c4762a1bSJed Brown      -da_processors_y <MY> number of processors in x direction
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petscdm.h>
12c4762a1bSJed Brown #include <petscdmda.h>
13c4762a1bSJed Brown 
FillLocalSubdomain(DM da,Vec gvec)14d71ae5a4SJacob Faibussowitsch PetscErrorCode FillLocalSubdomain(DM da, Vec gvec)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown   DMDALocalInfo info;
17c4762a1bSJed Brown   PetscMPIInt   rank;
18c4762a1bSJed Brown   PetscInt      i, j, k, l;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
229566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   if (info.dim == 3) {
25c4762a1bSJed Brown     PetscScalar ***g;
269566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, gvec, &g));
27c4762a1bSJed Brown     /* loop over ghosts */
28c4762a1bSJed Brown     for (k = info.zs; k < info.zs + info.zm; k++) {
29c4762a1bSJed Brown       for (j = info.ys; j < info.ys + info.ym; j++) {
30c4762a1bSJed Brown         for (i = info.xs; i < info.xs + info.xm; i++) {
31c4762a1bSJed Brown           g[k][j][info.dof * i + 0] = i;
32c4762a1bSJed Brown           g[k][j][info.dof * i + 1] = j;
33c4762a1bSJed Brown           g[k][j][info.dof * i + 2] = k;
34c4762a1bSJed Brown         }
35c4762a1bSJed Brown       }
36c4762a1bSJed Brown     }
379566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, gvec, &g));
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   if (info.dim == 2) {
40c4762a1bSJed Brown     PetscScalar **g;
419566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, gvec, &g));
42c4762a1bSJed Brown     /* loop over ghosts */
43c4762a1bSJed Brown     for (j = info.ys; j < info.ys + info.ym; j++) {
44c4762a1bSJed Brown       for (i = info.xs; i < info.xs + info.xm; i++) {
45c4762a1bSJed Brown         for (l = 0; l < info.dof; l++) {
46c4762a1bSJed Brown           g[j][info.dof * i + 0] = i;
47c4762a1bSJed Brown           g[j][info.dof * i + 1] = j;
48c4762a1bSJed Brown           g[j][info.dof * i + 2] = rank;
49c4762a1bSJed Brown         }
50c4762a1bSJed Brown       }
51c4762a1bSJed Brown     }
529566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, gvec, &g));
53c4762a1bSJed Brown   }
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
main(int argc,char ** argv)57d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
58d71ae5a4SJacob Faibussowitsch {
59c4762a1bSJed Brown   DM            da, *subda;
6072d683eeSStefano Zampini   PetscInt      i, dim = 3, dof = 3;
613b5e53cdSSajid Ali   PetscInt      M = 25, N = 25, P = 25;
62c4762a1bSJed Brown   PetscMPIInt   size, rank;
63c4762a1bSJed Brown   Vec           v;
64c4762a1bSJed Brown   Vec           slvec, sgvec;
65c4762a1bSJed Brown   IS           *ois, *iis;
66c4762a1bSJed Brown   VecScatter    oscata;
67c4762a1bSJed Brown   VecScatter   *iscat, *oscat, *gscat;
68c4762a1bSJed Brown   DMDALocalInfo info;
693b5e53cdSSajid Ali   PetscBool     patchis_offproc = PETSC_TRUE;
70c4762a1bSJed Brown 
71327415f7SBarry Smith   PetscFunctionBeginUser;
72*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
7472d683eeSStefano Zampini   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dof", &dof, 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) {
8072d683eeSStefano Zampini     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, 1, NULL, NULL, &da));
81c4762a1bSJed Brown   } else if (dim == 3) {
8272d683eeSStefano Zampini     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, dof, 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 
9172d683eeSStefano Zampini   /* TODO: broken for dim=2 */
9272d683eeSStefano Zampini   if (dim == 3) {
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 
1009566063dSJacob Faibussowitsch     PetscCall(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 */
1109566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(da, &lower, &upper, &patchis, patchis_offproc));
1119566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(da, &largevec));
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &smallvec));
1149566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(smallvec, info.dof * (upper.i - lower.i) * (upper.j - lower.j) * (upper.k - lower.k), PETSC_DECIDE));
1159566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(smallvec));
1169566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(smallvec, NULL, largevec, patchis, &patchscat));
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch     PetscCall(FillLocalSubdomain(subda[0], smallvec));
1199566063dSJacob Faibussowitsch     PetscCall(VecSet(largevec, 0));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(patchscat, smallvec, largevec, ADD_VALUES, SCATTER_FORWARD));
1229566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(patchscat, smallvec, largevec, ADD_VALUES, SCATTER_FORWARD));
1239566063dSJacob Faibussowitsch     PetscCall(ISView(patchis, PETSC_VIEWER_STDOUT_WORLD));
1249566063dSJacob Faibussowitsch     PetscCall(VecScatterView(patchscat, PETSC_VIEWER_STDOUT_WORLD));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown     for (i = 0; i < size; i++) {
12748a46eb9SPierre Jolivet       if (i == rank) PetscCall(VecView(smallvec, PETSC_VIEWER_STDOUT_SELF));
1289566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
129c4762a1bSJed Brown     }
130c4762a1bSJed Brown 
1319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
1329566063dSJacob Faibussowitsch     PetscCall(VecView(largevec, PETSC_VIEWER_STDOUT_WORLD));
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&smallvec));
1359566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(da, &largevec));
1369566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&patchis));
1379566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&patchscat));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* view the various parts */
141c4762a1bSJed Brown   {
142c4762a1bSJed Brown     for (i = 0; i < size; i++) {
143c4762a1bSJed Brown       if (i == rank) {
1449566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor %d: \n", i));
1459566063dSJacob Faibussowitsch         PetscCall(DMView(subda[0], PETSC_VIEWER_STDOUT_SELF));
146c4762a1bSJed Brown       }
1479566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(subda[0], &slvec));
1519566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(subda[0], &sgvec));
1529566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(da, &v));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown     /* test filling outer between the big DM and the small ones with the IS scatter*/
1559566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(v, ois[0], sgvec, NULL, &oscata));
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch     PetscCall(FillLocalSubdomain(subda[0], sgvec));
158c4762a1bSJed Brown 
1599566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscata, sgvec, v, ADD_VALUES, SCATTER_REVERSE));
1609566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscata, sgvec, v, ADD_VALUES, SCATTER_REVERSE));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown     /* test the local-to-local scatter */
163c4762a1bSJed Brown 
164c4762a1bSJed Brown     /* fill up the local subdomain and then add them together */
1659566063dSJacob Faibussowitsch     PetscCall(FillLocalSubdomain(da, v));
166c4762a1bSJed Brown 
1679566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat[0], v, slvec, ADD_VALUES, SCATTER_FORWARD));
1689566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat[0], v, slvec, ADD_VALUES, SCATTER_FORWARD));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown     /* test ghost scattering backwards */
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(VecSet(v, 0));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat[0], slvec, v, ADD_VALUES, SCATTER_REVERSE));
1779566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat[0], slvec, v, ADD_VALUES, SCATTER_REVERSE));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown     /* test overlap scattering backwards */
182c4762a1bSJed Brown 
1839566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(subda[0], slvec, ADD_VALUES, sgvec));
1849566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(subda[0], slvec, ADD_VALUES, sgvec));
185c4762a1bSJed Brown 
1869566063dSJacob Faibussowitsch     PetscCall(VecSet(v, 0));
187c4762a1bSJed Brown 
1889566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE));
1899566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE));
190c4762a1bSJed Brown 
1919566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
192c4762a1bSJed Brown 
193c4762a1bSJed Brown     /* test interior scattering backwards */
194c4762a1bSJed Brown 
1959566063dSJacob Faibussowitsch     PetscCall(VecSet(v, 0));
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(iscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE));
1989566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(iscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE));
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown     /* test matrix allocation */
203c4762a1bSJed Brown     for (i = 0; i < size; i++) {
204c4762a1bSJed Brown       if (i == rank) {
205c4762a1bSJed Brown         Mat m;
2069566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor %d: \n", i));
2079566063dSJacob Faibussowitsch         PetscCall(DMSetMatType(subda[0], MATAIJ));
2089566063dSJacob Faibussowitsch         PetscCall(DMCreateMatrix(subda[0], &m));
2099566063dSJacob Faibussowitsch         PetscCall(MatView(m, PETSC_VIEWER_STDOUT_SELF));
2109566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&m));
211c4762a1bSJed Brown       }
2129566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
213c4762a1bSJed Brown     }
2149566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(subda[0], &slvec));
2159566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(subda[0], &sgvec));
2169566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(da, &v));
217c4762a1bSJed Brown   }
218c4762a1bSJed Brown 
2199566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&subda[0]));
2209566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ois[0]));
2219566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iis[0]));
222c4762a1bSJed Brown 
2239566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&iscat[0]));
2249566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&oscat[0]));
2259566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&gscat[0]));
2269566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&oscata));
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(PetscFree(iscat));
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(oscat));
2309566063dSJacob Faibussowitsch   PetscCall(PetscFree(gscat));
2319566063dSJacob Faibussowitsch   PetscCall(PetscFree(oscata));
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree(subda));
2349566063dSJacob Faibussowitsch   PetscCall(PetscFree(ois));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFree(iis));
236c4762a1bSJed Brown 
2379566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2389566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
239b122ec5aSJacob Faibussowitsch   return 0;
240c4762a1bSJed Brown }
241