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 159371c9d4SSatish Balay PetscErrorCode FillLocalSubdomain(DM da, Vec gvec) { 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 } 54c4762a1bSJed Brown PetscFunctionReturn(0); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 579371c9d4SSatish Balay int main(int argc, char **argv) { 58c4762a1bSJed Brown DM da, *subda; 59c4762a1bSJed Brown PetscInt i, dim = 3; 603b5e53cdSSajid Ali PetscInt M = 25, N = 25, P = 25; 61c4762a1bSJed Brown PetscMPIInt size, rank; 62c4762a1bSJed Brown Vec v; 63c4762a1bSJed Brown Vec slvec, sgvec; 64c4762a1bSJed Brown IS *ois, *iis; 65c4762a1bSJed Brown VecScatter oscata; 66c4762a1bSJed Brown VecScatter *iscat, *oscat, *gscat; 67c4762a1bSJed Brown DMDALocalInfo info; 683b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 69c4762a1bSJed Brown 70327415f7SBarry Smith PetscFunctionBeginUser; 719566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Create distributed array and get vectors */ 759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 77c4762a1bSJed Brown if (dim == 2) { 789566063dSJacob 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)); 79c4762a1bSJed Brown } else if (dim == 3) { 809566063dSJacob 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)); 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 839566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 849566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(da, NULL, NULL, &iis, &ois, &subda)); 879566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecompositionScatters(da, 1, subda, &iscat, &oscat, &gscat)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown { 90c4762a1bSJed Brown DMDALocalInfo subinfo; 91c4762a1bSJed Brown MatStencil lower, upper; 923b5e53cdSSajid Ali IS patchis; 93c4762a1bSJed Brown Vec smallvec; 94c4762a1bSJed Brown Vec largevec; 95c4762a1bSJed Brown VecScatter patchscat; 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subda[0], &subinfo)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown lower.i = info.xs; 100c4762a1bSJed Brown lower.j = info.ys; 101c4762a1bSJed Brown lower.k = info.zs; 102c4762a1bSJed Brown upper.i = info.xs + info.xm; 103c4762a1bSJed Brown upper.j = info.ys + info.ym; 104c4762a1bSJed Brown upper.k = info.zs + info.zm; 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* test the patch IS as a thing to scatter to/from */ 1079566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da, &lower, &upper, &patchis, patchis_offproc)); 1089566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &largevec)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &smallvec)); 1119566063dSJacob Faibussowitsch PetscCall(VecSetSizes(smallvec, info.dof * (upper.i - lower.i) * (upper.j - lower.j) * (upper.k - lower.k), PETSC_DECIDE)); 1129566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(smallvec)); 1139566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(smallvec, NULL, largevec, patchis, &patchscat)); 114c4762a1bSJed Brown 1159566063dSJacob Faibussowitsch PetscCall(FillLocalSubdomain(subda[0], smallvec)); 1169566063dSJacob Faibussowitsch PetscCall(VecSet(largevec, 0)); 117c4762a1bSJed Brown 1189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(patchscat, smallvec, largevec, ADD_VALUES, SCATTER_FORWARD)); 1199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(patchscat, smallvec, largevec, ADD_VALUES, SCATTER_FORWARD)); 1209566063dSJacob Faibussowitsch PetscCall(ISView(patchis, PETSC_VIEWER_STDOUT_WORLD)); 1219566063dSJacob Faibussowitsch PetscCall(VecScatterView(patchscat, PETSC_VIEWER_STDOUT_WORLD)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown for (i = 0; i < size; i++) { 124*48a46eb9SPierre Jolivet if (i == rank) PetscCall(VecView(smallvec, PETSC_VIEWER_STDOUT_SELF)); 1259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 1299566063dSJacob Faibussowitsch PetscCall(VecView(largevec, PETSC_VIEWER_STDOUT_WORLD)); 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&smallvec)); 1329566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &largevec)); 1339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&patchis)); 1349566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&patchscat)); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* view the various parts */ 138c4762a1bSJed Brown { 139c4762a1bSJed Brown for (i = 0; i < size; i++) { 140c4762a1bSJed Brown if (i == rank) { 1419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor %d: \n", i)); 1429566063dSJacob Faibussowitsch PetscCall(DMView(subda[0], PETSC_VIEWER_STDOUT_SELF)); 143c4762a1bSJed Brown } 1449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(subda[0], &slvec)); 1489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(subda[0], &sgvec)); 1499566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &v)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* test filling outer between the big DM and the small ones with the IS scatter*/ 1529566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(v, ois[0], sgvec, NULL, &oscata)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(FillLocalSubdomain(subda[0], sgvec)); 155c4762a1bSJed Brown 1569566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscata, sgvec, v, ADD_VALUES, SCATTER_REVERSE)); 1579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscata, sgvec, v, ADD_VALUES, SCATTER_REVERSE)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* test the local-to-local scatter */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* fill up the local subdomain and then add them together */ 1629566063dSJacob Faibussowitsch PetscCall(FillLocalSubdomain(da, v)); 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat[0], v, slvec, ADD_VALUES, SCATTER_FORWARD)); 1659566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat[0], v, slvec, ADD_VALUES, SCATTER_FORWARD)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* test ghost scattering backwards */ 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(VecSet(v, 0)); 172c4762a1bSJed Brown 1739566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat[0], slvec, v, ADD_VALUES, SCATTER_REVERSE)); 1749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat[0], slvec, v, ADD_VALUES, SCATTER_REVERSE)); 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* test overlap scattering backwards */ 179c4762a1bSJed Brown 1809566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(subda[0], slvec, ADD_VALUES, sgvec)); 1819566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(subda[0], slvec, ADD_VALUES, sgvec)); 182c4762a1bSJed Brown 1839566063dSJacob Faibussowitsch PetscCall(VecSet(v, 0)); 184c4762a1bSJed Brown 1859566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE)); 1869566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE)); 187c4762a1bSJed Brown 1889566063dSJacob Faibussowitsch PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* test interior scattering backwards */ 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(VecSet(v, 0)); 193c4762a1bSJed Brown 1949566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(iscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE)); 1959566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(iscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE)); 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD)); 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* test matrix allocation */ 200c4762a1bSJed Brown for (i = 0; i < size; i++) { 201c4762a1bSJed Brown if (i == rank) { 202c4762a1bSJed Brown Mat m; 2039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor %d: \n", i)); 2049566063dSJacob Faibussowitsch PetscCall(DMSetMatType(subda[0], MATAIJ)); 2059566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(subda[0], &m)); 2069566063dSJacob Faibussowitsch PetscCall(MatView(m, PETSC_VIEWER_STDOUT_SELF)); 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&m)); 208c4762a1bSJed Brown } 2099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 210c4762a1bSJed Brown } 2119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(subda[0], &slvec)); 2129566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(subda[0], &sgvec)); 2139566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &v)); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 2169566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subda[0])); 2179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ois[0])); 2189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iis[0])); 219c4762a1bSJed Brown 2209566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&iscat[0])); 2219566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&oscat[0])); 2229566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&gscat[0])); 2239566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&oscata)); 224c4762a1bSJed Brown 2259566063dSJacob Faibussowitsch PetscCall(PetscFree(iscat)); 2269566063dSJacob Faibussowitsch PetscCall(PetscFree(oscat)); 2279566063dSJacob Faibussowitsch PetscCall(PetscFree(gscat)); 2289566063dSJacob Faibussowitsch PetscCall(PetscFree(oscata)); 229c4762a1bSJed Brown 2309566063dSJacob Faibussowitsch PetscCall(PetscFree(subda)); 2319566063dSJacob Faibussowitsch PetscCall(PetscFree(ois)); 2329566063dSJacob Faibussowitsch PetscCall(PetscFree(iis)); 233c4762a1bSJed Brown 2349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 2359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 236b122ec5aSJacob Faibussowitsch return 0; 237c4762a1bSJed Brown } 238