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\
5da81f932SPierre Jolivet sliceid - set the location where the slice will be extracted 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
main(int argc,char ** argv)15d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
16d71ae5a4SJacob Faibussowitsch {
173b5e53cdSSajid Ali PetscMPIInt rank, size; /* MPI rank and size */
183b5e53cdSSajid Ali PetscInt mx = 4, my = 4, mz = 4; /* Dimensions of parent vector */
193b5e53cdSSajid Ali PetscInt sliceid = 2; /* k (z) index to pick the slice */
203b5e53cdSSajid Ali PetscInt sliceaxis = 2; /* Select axis along which the slice will be extracted */
213b5e53cdSSajid Ali PetscInt dim = 3; /* Dimension of the parent vector */
223b5e53cdSSajid Ali PetscInt i, j, k; /* Iteration indices */
233b5e53cdSSajid Ali PetscInt ixs, iys, izs; /* Corner indices for 3D vector */
243b5e53cdSSajid Ali PetscInt ixm, iym, izm; /* Widths of parent vector */
253b5e53cdSSajid Ali PetscScalar ***vecdata3d; /* Pointer to access 3d parent vector */
263b5e53cdSSajid Ali PetscScalar **vecdata2d; /* Pointer to access 2d parent vector */
273b5e53cdSSajid Ali DM da; /* 2D/3D DMDA object */
283b5e53cdSSajid Ali Vec vec_full; /* Parent vector */
293b5e53cdSSajid Ali Vec vec_slice; /* Slice vector */
303b5e53cdSSajid Ali MatStencil lower, upper; /* Stencils to select slice */
313b5e53cdSSajid Ali IS selectis; /* IS to select slice and extract subvector */
323b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_FALSE; /* flag to DMDACreatePatchIS indicating that off-proc values are to be ignored */
333b5e53cdSSajid Ali
343b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
353b5e53cdSSajid Ali Initialize program and set problem parameters
363b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37327415f7SBarry Smith PetscFunctionBeginUser;
38*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
413b5e53cdSSajid Ali
429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, NULL));
439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, NULL));
449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, NULL));
459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-sliceid", &sliceid, NULL));
479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-sliceaxis", &sliceaxis, NULL));
483b5e53cdSSajid Ali
493b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
503b5e53cdSSajid Ali Create DMDA object.
513b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
523b5e53cdSSajid Ali if (dim == 3) {
539371c9d4SSatish 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));
549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
559566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
563b5e53cdSSajid Ali } else {
579371c9d4SSatish 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));
589566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
599566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
603b5e53cdSSajid Ali }
613b5e53cdSSajid Ali
623b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
633b5e53cdSSajid Ali Create the parent vector
643b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
659566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &vec_full));
669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec_full, "full_vector"));
673b5e53cdSSajid Ali
683b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
693b5e53cdSSajid Ali Populate the 3D vector
703b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
719566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &ixs, &iys, &izs, &ixm, &iym, &izm));
723b5e53cdSSajid Ali if (dim == 3) {
739566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, vec_full, &vecdata3d));
743b5e53cdSSajid Ali for (k = izs; k < izs + izm; k++) {
753b5e53cdSSajid Ali for (j = iys; j < iys + iym; j++) {
76ad540459SPierre Jolivet for (i = ixs; i < ixs + ixm; i++) vecdata3d[k][j][i] = ((i - mx / 2) * (j + mx / 2)) + k * 100;
773b5e53cdSSajid Ali }
783b5e53cdSSajid Ali }
799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, vec_full, &vecdata3d));
803b5e53cdSSajid Ali } else {
819566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, vec_full, &vecdata2d));
823b5e53cdSSajid Ali for (j = iys; j < iys + iym; j++) {
83ad540459SPierre Jolivet for (i = ixs; i < ixs + ixm; i++) vecdata2d[j][i] = ((i - mx / 2) * (j + mx / 2));
843b5e53cdSSajid Ali }
859566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, vec_full, &vecdata2d));
863b5e53cdSSajid Ali }
873b5e53cdSSajid Ali
883b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
893b5e53cdSSajid Ali Get an IS corresponding to a 2D slice
903b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
913b5e53cdSSajid Ali if (sliceaxis == 0) {
929371c9d4SSatish Balay lower.i = sliceid;
939371c9d4SSatish Balay lower.j = 0;
949371c9d4SSatish Balay lower.k = 0;
959371c9d4SSatish Balay lower.c = 1;
969371c9d4SSatish Balay upper.i = sliceid;
979371c9d4SSatish Balay upper.j = my;
989371c9d4SSatish Balay upper.k = mz;
999371c9d4SSatish Balay upper.c = 1;
1003b5e53cdSSajid Ali } else if (sliceaxis == 1) {
1019371c9d4SSatish Balay lower.i = 0;
1029371c9d4SSatish Balay lower.j = sliceid;
1039371c9d4SSatish Balay lower.k = 0;
1049371c9d4SSatish Balay lower.c = 1;
1059371c9d4SSatish Balay upper.i = mx;
1069371c9d4SSatish Balay upper.j = sliceid;
1079371c9d4SSatish Balay upper.k = mz;
1089371c9d4SSatish Balay upper.c = 1;
1093b5e53cdSSajid Ali } else if (sliceaxis == 2 && dim == 3) {
1109371c9d4SSatish Balay lower.i = 0;
1119371c9d4SSatish Balay lower.j = 0;
1129371c9d4SSatish Balay lower.k = sliceid;
1139371c9d4SSatish Balay lower.c = 1;
1149371c9d4SSatish Balay upper.i = mx;
1159371c9d4SSatish Balay upper.j = my;
1169371c9d4SSatish Balay upper.k = sliceid;
1179371c9d4SSatish Balay upper.c = 1;
1183b5e53cdSSajid Ali }
1199566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da, &lower, &upper, &selectis, patchis_offproc));
1209566063dSJacob Faibussowitsch PetscCall(ISView(selectis, PETSC_VIEWER_STDOUT_WORLD));
1213b5e53cdSSajid Ali
1223b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1233b5e53cdSSajid Ali Use the obtained IS to extract the slice as a subvector
1243b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1259566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vec_full, selectis, &vec_slice));
1263b5e53cdSSajid Ali
1273b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1283b5e53cdSSajid Ali View the extracted subvector
1293b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1309566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_DENSE));
1319566063dSJacob Faibussowitsch PetscCall(VecView(vec_slice, PETSC_VIEWER_STDOUT_WORLD));
1329566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
1333b5e53cdSSajid Ali
1343b5e53cdSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1353b5e53cdSSajid Ali Restore subvector, destroy data structures and exit.
1363b5e53cdSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1379566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vec_full, selectis, &vec_slice));
1383b5e53cdSSajid Ali
1399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&selectis));
1409566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_full));
1423b5e53cdSSajid Ali
1439566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
144b122ec5aSJacob Faibussowitsch return 0;
1453b5e53cdSSajid Ali }
1463b5e53cdSSajid Ali
1473b5e53cdSSajid Ali /*TEST
1483b5e53cdSSajid Ali
1493b5e53cdSSajid Ali test:
1503b5e53cdSSajid Ali nsize: 1
1513b5e53cdSSajid Ali args: -sliceaxis 0
1523b5e53cdSSajid Ali
1533b5e53cdSSajid Ali test:
1543b5e53cdSSajid Ali suffix: 2
1553b5e53cdSSajid Ali nsize: 2
1563b5e53cdSSajid Ali args: -sliceaxis 1
1573b5e53cdSSajid Ali
1583b5e53cdSSajid Ali test:
1593b5e53cdSSajid Ali suffix: 3
1603b5e53cdSSajid Ali nsize: 4
1613b5e53cdSSajid Ali args: -sliceaxis 2
1623b5e53cdSSajid Ali
1633b5e53cdSSajid Ali test:
1643b5e53cdSSajid Ali suffix: 4
1653b5e53cdSSajid Ali nsize: 2
1663b5e53cdSSajid Ali args: -sliceaxis 1 -dim 2
1673b5e53cdSSajid Ali
1683b5e53cdSSajid Ali test:
1693b5e53cdSSajid Ali suffix: 5
1703b5e53cdSSajid Ali nsize: 3
1713b5e53cdSSajid Ali args: -sliceaxis 0 -dim 2
1723b5e53cdSSajid Ali
1733b5e53cdSSajid Ali TEST*/
174