187adc40fSSajid Ali static char help[] = "Extract a 2D slice in natural ordering from a 3D vector, Command line options : \n\
287adc40fSSajid Ali Mx/My/Mz - set the dimensions of the parent vector \n\
387adc40fSSajid Ali sliceaxis - string describing the axis along which the sice will be selected : X, Y, Z \n\
487adc40fSSajid Ali gp - global grid point number along the sliceaxis direction where the slice will be extracted from the parent vector \n";
587adc40fSSajid Ali
687adc40fSSajid Ali /*
787adc40fSSajid Ali This example shows to extract a 2D slice in natural ordering
887adc40fSSajid Ali from a 3D DMDA vector (first by extracting the slice and then
987adc40fSSajid Ali by converting it to natural ordering)
1087adc40fSSajid Ali */
1187adc40fSSajid Ali
1287adc40fSSajid Ali #include <petscdmda.h>
1387adc40fSSajid Ali
1487adc40fSSajid Ali const char *const sliceaxes[] = {"X", "Y", "Z", "sliceaxis", "DM_", NULL};
1587adc40fSSajid Ali
main(int argc,char ** argv)16d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
17d71ae5a4SJacob Faibussowitsch {
1887adc40fSSajid Ali DM da3D; /* 3D DMDA object */
1987adc40fSSajid Ali DM da2D; /* 2D DMDA object */
2087adc40fSSajid Ali Vec vec_full; /* Parent vector */
2187adc40fSSajid Ali Vec vec_extracted; /* Extracted slice vector (in DMDA ordering) */
2287adc40fSSajid Ali Vec vec_slice; /* vector in natural ordering */
2387adc40fSSajid Ali Vec vec_slice_g; /* aliased vector in natural ordering */
2487adc40fSSajid Ali IS patchis_3d; /* IS to select slice and extract subvector */
2587adc40fSSajid Ali IS patchis_2d; /* Patch IS for 2D vector, will be converted to application ordering */
2687adc40fSSajid Ali IS scatis_extracted_slice; /* PETSc indexed IS for extracted slice */
2787adc40fSSajid Ali IS scatis_natural_slice; /* natural/application ordered IS for slice*/
2887adc40fSSajid Ali IS scatis_natural_slice_g; /* aliased natural/application ordered IS for slice */
2987adc40fSSajid Ali VecScatter vscat; /* scatter slice in DMDA ordering <-> slice in column major ordering */
3087adc40fSSajid Ali AO da2D_ao; /* AO associated with 2D DMDA */
3187adc40fSSajid Ali MPI_Comm subset_mpi_comm = MPI_COMM_NULL; /* MPI communicator where the slice lives */
3287adc40fSSajid Ali PetscScalar ***vecdata3d; /* Pointer to access 3d parent vector */
3387adc40fSSajid Ali const PetscScalar *array; /* pointer to create aliased Vec */
3487adc40fSSajid Ali PetscInt Mx = 4, My = 4, Mz = 4; /* Dimensions for 3D DMDA */
3587adc40fSSajid Ali const PetscInt *l1, *l2; /* 3D DMDA layout */
3687adc40fSSajid Ali PetscInt M1 = -1, M2 = -1; /* Dimensions for 2D DMDA */
3787adc40fSSajid Ali PetscInt m1 = -1, m2 = -1; /* Layouts for 2D DMDA */
3887adc40fSSajid Ali PetscInt gp = 2; /* grid point along sliceaxis to pick the slice */
3987adc40fSSajid Ali DMDirection sliceaxis = DM_X; /* Select axis along which the slice will be extracted */
4087adc40fSSajid Ali PetscInt i, j, k; /* Iteration indices */
4187adc40fSSajid Ali PetscInt ixs, iys, izs; /* Corner indices for 3D vector */
4287adc40fSSajid Ali PetscInt ixm, iym, izm; /* Widths of parent vector */
4387adc40fSSajid Ali PetscInt low, high; /* ownership range indices */
4487adc40fSSajid Ali PetscInt in; /* local size index for IS*/
4587adc40fSSajid Ali PetscInt vn; /* local size index */
4687adc40fSSajid Ali const PetscInt *is_array; /* pointer to create aliased IS */
4787adc40fSSajid Ali MatStencil lower, upper; /* Stencils to select slice for Vec */
4887adc40fSSajid Ali PetscBool patchis_offproc = PETSC_FALSE; /* flag to DMDACreatePatchIS indicating that off-proc values are to be ignored */
4987adc40fSSajid Ali PetscMPIInt rank, size; /* MPI rank and size */
5087adc40fSSajid Ali
5187adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5287adc40fSSajid Ali Initialize program and set problem parameters
5387adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54327415f7SBarry Smith PetscFunctionBeginUser;
55*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
5887adc40fSSajid Ali
59d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "ex22 DMDA tutorial example options", "DMDA");
601690c2aeSBarry Smith PetscCall(PetscOptionsRangeInt("-Mx", "dimension along x-axis", "ex22.c", Mx, &Mx, NULL, 0, PETSC_INT_MAX));
611690c2aeSBarry Smith PetscCall(PetscOptionsRangeInt("-My", "dimension along y-axis", "ex22.c", My, &My, NULL, 0, PETSC_INT_MAX));
621690c2aeSBarry Smith PetscCall(PetscOptionsRangeInt("-Mz", "dimension along z-axis", "ex22.c", Mz, &Mz, NULL, 0, PETSC_INT_MAX));
639566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sliceaxis", "axis along which 2D slice is extracted from : X, Y, Z", "", sliceaxes, (PetscEnum)sliceaxis, (PetscEnum *)&sliceaxis, NULL));
641690c2aeSBarry Smith PetscCall(PetscOptionsRangeInt("-gp", "index along sliceaxis at which 2D slice is extracted", "ex22.c", gp, &gp, NULL, 0, PETSC_INT_MAX));
65d0609cedSBarry Smith PetscOptionsEnd();
6687adc40fSSajid Ali
6787adc40fSSajid Ali /* Ensure that the requested slice is not out of bounds for the selected axis */
6887adc40fSSajid Ali if (sliceaxis == DM_X) {
6908401ef6SPierre Jolivet PetscCheck(gp <= Mx, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "grid point along sliceaxis is larger than largest index!");
7087adc40fSSajid Ali } else if (sliceaxis == DM_Y) {
7108401ef6SPierre Jolivet PetscCheck(gp <= My, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "grid point along sliceaxis is larger than largest index!");
7287adc40fSSajid Ali } else if (sliceaxis == DM_Z) {
7308401ef6SPierre Jolivet PetscCheck(gp <= Mz, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "grid point along sliceaxis is larger than largest index!");
7487adc40fSSajid Ali }
7587adc40fSSajid Ali
7687adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7787adc40fSSajid Ali Create 3D DMDA object.
7887adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
799371c9d4SSatish 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, &da3D));
809566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da3D));
819566063dSJacob Faibussowitsch PetscCall(DMSetUp(da3D));
8287adc40fSSajid Ali
8387adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8487adc40fSSajid Ali Create the parent vector
8587adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
869566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da3D, &vec_full));
879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec_full, "full_vector"));
8887adc40fSSajid Ali
8987adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9087adc40fSSajid Ali Populate the 3D vector
9187adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
929566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da3D, &ixs, &iys, &izs, &ixm, &iym, &izm));
939566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3D, vec_full, &vecdata3d));
9487adc40fSSajid Ali for (k = izs; k < izs + izm; k++) {
9587adc40fSSajid Ali for (j = iys; j < iys + iym; j++) {
96ad540459SPierre Jolivet for (i = ixs; i < ixs + ixm; i++) vecdata3d[k][j][i] = ((i - Mx / 2.0) * (j + Mx / 2.0)) + k * 100;
9787adc40fSSajid Ali }
9887adc40fSSajid Ali }
999566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3D, vec_full, &vecdata3d));
10087adc40fSSajid Ali
10187adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
10287adc40fSSajid Ali Get an IS corresponding to a 2D slice
10387adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
10487adc40fSSajid Ali if (sliceaxis == DM_X) {
1059371c9d4SSatish Balay lower.i = gp;
1069371c9d4SSatish Balay lower.j = 0;
1079371c9d4SSatish Balay lower.k = 0;
1089371c9d4SSatish Balay upper.i = gp;
1099371c9d4SSatish Balay upper.j = My;
1109371c9d4SSatish Balay upper.k = Mz;
11187adc40fSSajid Ali } else if (sliceaxis == DM_Y) {
1129371c9d4SSatish Balay lower.i = 0;
1139371c9d4SSatish Balay lower.j = gp;
1149371c9d4SSatish Balay lower.k = 0;
1159371c9d4SSatish Balay upper.i = Mx;
1169371c9d4SSatish Balay upper.j = gp;
1179371c9d4SSatish Balay upper.k = Mz;
11887adc40fSSajid Ali } else if (sliceaxis == DM_Z) {
1199371c9d4SSatish Balay lower.i = 0;
1209371c9d4SSatish Balay lower.j = 0;
1219371c9d4SSatish Balay lower.k = gp;
1229371c9d4SSatish Balay upper.i = Mx;
1239371c9d4SSatish Balay upper.j = My;
1249371c9d4SSatish Balay upper.k = gp;
12587adc40fSSajid Ali }
1269566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da3D, &lower, &upper, &patchis_3d, patchis_offproc));
1279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n IS to select slice from 3D DMDA vector : \n"));
1289566063dSJacob Faibussowitsch PetscCall(ISView(patchis_3d, PETSC_VIEWER_STDOUT_WORLD));
12987adc40fSSajid Ali
13087adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
13187adc40fSSajid Ali Use the obtained IS to extract the slice as a subvector
13287adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1339566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vec_full, patchis_3d, &vec_extracted));
13487adc40fSSajid Ali
13587adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
13687adc40fSSajid Ali View the extracted subvector
13787adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1389566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_DENSE));
1399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in DMDA ordering : \n"));
1409566063dSJacob Faibussowitsch PetscCall(VecView(vec_extracted, PETSC_VIEWER_STDOUT_WORLD));
1419566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
14287adc40fSSajid Ali
14387adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14487adc40fSSajid Ali Query 3D DMDA layout, get the subset MPI communicator
14587adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
14687adc40fSSajid Ali if (sliceaxis == DM_X) {
1479566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL));
1489566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da3D, NULL, &l1, &l2));
1499371c9d4SSatish Balay M1 = My;
1509371c9d4SSatish Balay M2 = Mz;
1519566063dSJacob Faibussowitsch PetscCall(DMDAGetProcessorSubset(da3D, DM_X, gp, &subset_mpi_comm));
15287adc40fSSajid Ali } else if (sliceaxis == DM_Y) {
1539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, NULL, &m2, NULL, NULL, NULL, NULL, NULL, NULL));
1549566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da3D, &l1, NULL, &l2));
1559371c9d4SSatish Balay M1 = Mx;
1569371c9d4SSatish Balay M2 = Mz;
1579566063dSJacob Faibussowitsch PetscCall(DMDAGetProcessorSubset(da3D, DM_Y, gp, &subset_mpi_comm));
15887adc40fSSajid Ali } else if (sliceaxis == DM_Z) {
1599566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
1609566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da3D, &l1, &l2, NULL));
1619371c9d4SSatish Balay M1 = Mx;
1629371c9d4SSatish Balay M2 = My;
1639566063dSJacob Faibussowitsch PetscCall(DMDAGetProcessorSubset(da3D, DM_Z, gp, &subset_mpi_comm));
16487adc40fSSajid Ali }
16587adc40fSSajid Ali
16687adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
16787adc40fSSajid Ali Create 2D DMDA object,
16887adc40fSSajid Ali vector (that will hold the slice as a column major flattened array) &
16987adc40fSSajid Ali index set (that will be used for scattering to the column major
17087adc40fSSajid Ali indexed slice vector)
17187adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
17287adc40fSSajid Ali if (subset_mpi_comm != MPI_COMM_NULL) {
1739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subset_mpi_comm, &size));
1749566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(subset_mpi_comm, "subset MPI subcomm size is : %d, includes global rank : %d \n", size, rank));
1759566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(subset_mpi_comm, PETSC_STDOUT));
176d0609cedSBarry Smith PetscCall(DMDACreate2d(subset_mpi_comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M1, M2, m1, m2, 1, 1, l1, l2, &da2D));
1779566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da2D));
1789566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2D));
17987adc40fSSajid Ali
18087adc40fSSajid Ali /* Create a 2D patch IS for the slice */
1819371c9d4SSatish Balay lower.i = 0;
1829371c9d4SSatish Balay lower.j = 0;
1839371c9d4SSatish Balay upper.i = M1;
1849371c9d4SSatish Balay upper.j = M2;
1859566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da2D, &lower, &upper, &patchis_2d, patchis_offproc));
18687adc40fSSajid Ali
18787adc40fSSajid Ali /* Convert the 2D patch IS to natural indexing (column major flattened) */
1889566063dSJacob Faibussowitsch PetscCall(ISDuplicate(patchis_2d, &scatis_natural_slice));
1899566063dSJacob Faibussowitsch PetscCall(DMDAGetAO(da2D, &da2D_ao));
1909566063dSJacob Faibussowitsch PetscCall(AOPetscToApplicationIS(da2D_ao, scatis_natural_slice));
1919566063dSJacob Faibussowitsch PetscCall(ISGetIndices(scatis_natural_slice, &is_array));
1929566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(scatis_natural_slice, &in));
19387adc40fSSajid Ali
19487adc40fSSajid Ali /* Create an aliased IS on the 3D DMDA's communicator */
1959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, in, is_array, PETSC_USE_POINTER, &scatis_natural_slice_g));
1969566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(scatis_natural_slice, &is_array));
19787adc40fSSajid Ali
19887adc40fSSajid Ali /* Create a 2D DMDA global vector */
1999566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da2D, &vec_slice));
2009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec_slice, "slice_vector_natural"));
2019566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec_slice, &vn));
2029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vec_slice, &array));
20387adc40fSSajid Ali
20487adc40fSSajid Ali /* Create an aliased version of the above on the 3D DMDA's communicator */
2059566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, vn, M1 * M2, array, &vec_slice_g));
2069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vec_slice, &array));
20787adc40fSSajid Ali } else {
20887adc40fSSajid Ali /* Ranks not part of the subset MPI communicator provide no entries, but the routines for creating
20987adc40fSSajid Ali the IS and Vec on the 3D DMDA's communicator still need to called, since they are collective routines */
2109566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 0, NULL, PETSC_USE_POINTER, &scatis_natural_slice_g));
2119566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, 0, M1 * M2, NULL, &vec_slice_g));
21287adc40fSSajid Ali }
2139566063dSJacob Faibussowitsch PetscCall(PetscBarrier(NULL));
21487adc40fSSajid Ali
21587adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
21687adc40fSSajid Ali Create IS that maps from the extracted slice vector
21787adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2189566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vec_extracted, &low, &high));
2199566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, high - low, low, 1, &scatis_extracted_slice));
22087adc40fSSajid Ali
22187adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22287adc40fSSajid Ali Scatter extracted subvector -> natural 2D slice vector
22387adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2249566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec_extracted, scatis_extracted_slice, vec_slice_g, scatis_natural_slice_g, &vscat));
2259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD));
2269566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD));
22787adc40fSSajid Ali
22887adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22987adc40fSSajid Ali View the natural 2D slice vector
23087adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2319566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_DENSE));
2329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in natural ordering : \n"));
2339566063dSJacob Faibussowitsch PetscCall(VecView(vec_slice_g, PETSC_VIEWER_STDOUT_WORLD));
2349566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
23587adc40fSSajid Ali
23687adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23787adc40fSSajid Ali Restore subvector
23887adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2399566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vec_full, patchis_3d, &vec_extracted));
24087adc40fSSajid Ali
24187adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24287adc40fSSajid Ali Destroy data structures and exit.
24387adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_full));
2459566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat));
2469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&scatis_extracted_slice));
2479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&scatis_natural_slice_g));
2489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_slice_g));
2499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&patchis_3d));
2509566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da3D));
25187adc40fSSajid Ali
25287adc40fSSajid Ali if (subset_mpi_comm != MPI_COMM_NULL) {
2539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&patchis_2d));
2549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&scatis_natural_slice));
2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_slice));
2569566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da2D));
2579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subset_mpi_comm));
25887adc40fSSajid Ali }
25987adc40fSSajid Ali
2609566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
261b122ec5aSJacob Faibussowitsch return 0;
26287adc40fSSajid Ali }
26387adc40fSSajid Ali
26487adc40fSSajid Ali /*TEST
26587adc40fSSajid Ali
26687adc40fSSajid Ali test:
26787adc40fSSajid Ali nsize: 1
26887adc40fSSajid Ali args: -sliceaxis X -gp 0
26987adc40fSSajid Ali
27087adc40fSSajid Ali test:
27187adc40fSSajid Ali suffix: 2
27287adc40fSSajid Ali nsize: 2
27387adc40fSSajid Ali args: -sliceaxis Y -gp 1
2742071a92dSBarry Smith filter: grep -v "subset MPI subcomm"
27587adc40fSSajid Ali
27687adc40fSSajid Ali test:
27787adc40fSSajid Ali suffix: 3
27887adc40fSSajid Ali nsize: 3
27987adc40fSSajid Ali args: -sliceaxis Z -gp 2
2802071a92dSBarry Smith filter: grep -v "subset MPI subcomm"
28187adc40fSSajid Ali
28287adc40fSSajid Ali test:
28387adc40fSSajid Ali suffix: 4
28487adc40fSSajid Ali nsize: 4
28587adc40fSSajid Ali args: -sliceaxis X -gp 2
2862071a92dSBarry Smith filter: grep -v "subset MPI subcomm"
28787adc40fSSajid Ali
28887adc40fSSajid Ali test:
28987adc40fSSajid Ali suffix: 5
29087adc40fSSajid Ali nsize: 4
29187adc40fSSajid Ali args: -sliceaxis Z -gp 1
2922071a92dSBarry Smith filter: grep -v "subset MPI subcomm"
29387adc40fSSajid Ali
29487adc40fSSajid Ali TEST*/
295