13b5e53cdSSajid Ali static char help[] = "Use DMDACreatePatchIS to extract a slice from a vector, Command line options :\n\ 23b5e53cdSSajid Ali mx/my/mz - set the dimensions of the parent vector\n\ 33b5e53cdSSajid Ali dim - set the dimensionality of the parent vector (2,3)\n\ 43b5e53cdSSajid Ali sliceaxis - Integer describing the axis along which the sice will be selected (0-X, 1-Y, 2-Z)\n\ 53b5e53cdSSajid Ali sliceid - set the location where the slice will be extraced from the parent vector\n"; 63b5e53cdSSajid Ali 73b5e53cdSSajid Ali /* 83b5e53cdSSajid Ali This test checks the functionality of DMDACreatePatchIS when 93b5e53cdSSajid Ali extracting a 2D vector from a 3D vector and 1D vector from a 103b5e53cdSSajid Ali 2D vector. 113b5e53cdSSajid Ali */ 123b5e53cdSSajid Ali 133b5e53cdSSajid Ali #include <petscdmda.h> 143b5e53cdSSajid Ali 15*9371c9d4SSatish Balay int main(int argc, char **argv) { 163b5e53cdSSajid Ali PetscMPIInt rank, size; /* MPI rank and size */ 173b5e53cdSSajid Ali PetscInt mx = 4, my = 4, mz = 4; /* Dimensions of parent vector */ 183b5e53cdSSajid Ali PetscInt sliceid = 2; /* k (z) index to pick the slice */ 193b5e53cdSSajid Ali PetscInt sliceaxis = 2; /* Select axis along which the slice will be extracted */ 203b5e53cdSSajid Ali PetscInt dim = 3; /* Dimension of the parent vector */ 213b5e53cdSSajid Ali PetscInt i, j, k; /* Iteration indices */ 223b5e53cdSSajid Ali PetscInt ixs, iys, izs; /* Corner indices for 3D vector */ 233b5e53cdSSajid Ali PetscInt ixm, iym, izm; /* Widths of parent vector */ 243b5e53cdSSajid Ali PetscScalar ***vecdata3d; /* Pointer to access 3d parent vector */ 253b5e53cdSSajid Ali PetscScalar **vecdata2d; /* Pointer to access 2d parent vector */ 263b5e53cdSSajid Ali DM da; /* 2D/3D DMDA object */ 273b5e53cdSSajid Ali Vec vec_full; /* Parent vector */ 283b5e53cdSSajid Ali Vec vec_slice; /* Slice vector */ 293b5e53cdSSajid Ali MatStencil lower, upper; /* Stencils to select slice */ 303b5e53cdSSajid Ali IS selectis; /* IS to select slice and extract subvector */ 313b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_FALSE; /* flag to DMDACreatePatchIS indicating that off-proc values are to be ignored */ 323b5e53cdSSajid Ali 333b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 343b5e53cdSSajid Ali Initialize program and set problem parameters 353b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 36327415f7SBarry Smith PetscFunctionBeginUser; 379566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 403b5e53cdSSajid Ali 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, NULL)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-sliceid", &sliceid, NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-sliceaxis", &sliceaxis, NULL)); 473b5e53cdSSajid Ali 483b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 493b5e53cdSSajid Ali Create DMDA object. 503b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/ 513b5e53cdSSajid Ali if (dim == 3) { 52*9371c9d4SSatish Balay PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, mx, my, mz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &da)); 539566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 549566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 553b5e53cdSSajid Ali } else { 56*9371c9d4SSatish Balay PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, mx, my, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 589566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 593b5e53cdSSajid Ali } 603b5e53cdSSajid Ali 613b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 623b5e53cdSSajid Ali Create the parent vector 633b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/ 649566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &vec_full)); 659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec_full, "full_vector")); 663b5e53cdSSajid Ali 673b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 683b5e53cdSSajid Ali Populate the 3D vector 693b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 709566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &ixs, &iys, &izs, &ixm, &iym, &izm)); 713b5e53cdSSajid Ali if (dim == 3) { 729566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, vec_full, &vecdata3d)); 733b5e53cdSSajid Ali for (k = izs; k < izs + izm; k++) { 743b5e53cdSSajid Ali for (j = iys; j < iys + iym; j++) { 75*9371c9d4SSatish Balay for (i = ixs; i < ixs + ixm; i++) { vecdata3d[k][j][i] = ((i - mx / 2) * (j + mx / 2)) + k * 100; } 763b5e53cdSSajid Ali } 773b5e53cdSSajid Ali } 789566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, vec_full, &vecdata3d)); 793b5e53cdSSajid Ali } else { 809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, vec_full, &vecdata2d)); 813b5e53cdSSajid Ali for (j = iys; j < iys + iym; j++) { 82*9371c9d4SSatish Balay for (i = ixs; i < ixs + ixm; i++) { vecdata2d[j][i] = ((i - mx / 2) * (j + mx / 2)); } 833b5e53cdSSajid Ali } 849566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, vec_full, &vecdata2d)); 853b5e53cdSSajid Ali } 863b5e53cdSSajid Ali 873b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 883b5e53cdSSajid Ali Get an IS corresponding to a 2D slice 893b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 903b5e53cdSSajid Ali if (sliceaxis == 0) { 91*9371c9d4SSatish Balay lower.i = sliceid; 92*9371c9d4SSatish Balay lower.j = 0; 93*9371c9d4SSatish Balay lower.k = 0; 94*9371c9d4SSatish Balay lower.c = 1; 95*9371c9d4SSatish Balay upper.i = sliceid; 96*9371c9d4SSatish Balay upper.j = my; 97*9371c9d4SSatish Balay upper.k = mz; 98*9371c9d4SSatish Balay upper.c = 1; 993b5e53cdSSajid Ali } else if (sliceaxis == 1) { 100*9371c9d4SSatish Balay lower.i = 0; 101*9371c9d4SSatish Balay lower.j = sliceid; 102*9371c9d4SSatish Balay lower.k = 0; 103*9371c9d4SSatish Balay lower.c = 1; 104*9371c9d4SSatish Balay upper.i = mx; 105*9371c9d4SSatish Balay upper.j = sliceid; 106*9371c9d4SSatish Balay upper.k = mz; 107*9371c9d4SSatish Balay upper.c = 1; 1083b5e53cdSSajid Ali } else if (sliceaxis == 2 && dim == 3) { 109*9371c9d4SSatish Balay lower.i = 0; 110*9371c9d4SSatish Balay lower.j = 0; 111*9371c9d4SSatish Balay lower.k = sliceid; 112*9371c9d4SSatish Balay lower.c = 1; 113*9371c9d4SSatish Balay upper.i = mx; 114*9371c9d4SSatish Balay upper.j = my; 115*9371c9d4SSatish Balay upper.k = sliceid; 116*9371c9d4SSatish Balay upper.c = 1; 1173b5e53cdSSajid Ali } 1189566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da, &lower, &upper, &selectis, patchis_offproc)); 1199566063dSJacob Faibussowitsch PetscCall(ISView(selectis, PETSC_VIEWER_STDOUT_WORLD)); 1203b5e53cdSSajid Ali 1213b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1223b5e53cdSSajid Ali Use the obtained IS to extract the slice as a subvector 1233b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vec_full, selectis, &vec_slice)); 1253b5e53cdSSajid Ali 1263b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1273b5e53cdSSajid Ali View the extracted subvector 1283b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_DENSE)); 1309566063dSJacob Faibussowitsch PetscCall(VecView(vec_slice, PETSC_VIEWER_STDOUT_WORLD)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 1323b5e53cdSSajid Ali 1333b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1343b5e53cdSSajid Ali Restore subvector, destroy data structures and exit. 1353b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vec_full, selectis, &vec_slice)); 1373b5e53cdSSajid Ali 1389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selectis)); 1399566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_full)); 1413b5e53cdSSajid Ali 1429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 143b122ec5aSJacob Faibussowitsch return 0; 1443b5e53cdSSajid Ali } 1453b5e53cdSSajid Ali 1463b5e53cdSSajid Ali /*TEST 1473b5e53cdSSajid Ali 1483b5e53cdSSajid Ali test: 1493b5e53cdSSajid Ali nsize: 1 1503b5e53cdSSajid Ali args: -sliceaxis 0 1513b5e53cdSSajid Ali 1523b5e53cdSSajid Ali test: 1533b5e53cdSSajid Ali suffix: 2 1543b5e53cdSSajid Ali nsize: 2 1553b5e53cdSSajid Ali args: -sliceaxis 1 1563b5e53cdSSajid Ali 1573b5e53cdSSajid Ali test: 1583b5e53cdSSajid Ali suffix: 3 1593b5e53cdSSajid Ali nsize: 4 1603b5e53cdSSajid Ali args: -sliceaxis 2 1613b5e53cdSSajid Ali 1623b5e53cdSSajid Ali test: 1633b5e53cdSSajid Ali suffix: 4 1643b5e53cdSSajid Ali nsize: 2 1653b5e53cdSSajid Ali args: -sliceaxis 1 -dim 2 1663b5e53cdSSajid Ali 1673b5e53cdSSajid Ali test: 1683b5e53cdSSajid Ali suffix: 5 1693b5e53cdSSajid Ali nsize: 3 1703b5e53cdSSajid Ali args: -sliceaxis 0 -dim 2 1713b5e53cdSSajid Ali 1723b5e53cdSSajid Ali TEST*/ 173