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 1687adc40fSSajid Ali int main(int argc,char **argv) 1787adc40fSSajid Ali { 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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 54*327415f7SBarry Smith PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char*)0, 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"); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-Mx", "dimension along x-axis", "ex22.c", Mx, &Mx, NULL, 0, PETSC_MAX_INT)); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-My", "dimension along y-axis", "ex22.c", My, &My, NULL, 0, PETSC_MAX_INT)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-Mz", "dimension along z-axis", "ex22.c", Mz, &Mz, NULL, 0, PETSC_MAX_INT)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sliceaxis","axis along which 2D slice is extracted from : X, Y, Z","",sliceaxes,(PetscEnum)sliceaxis,(PetscEnum*)&sliceaxis,NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-gp", "index along sliceaxis at which 2D slice is extracted", "ex22.c", gp, &gp, NULL, 0, PETSC_MAX_INT)); 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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/ 79d0609cedSBarry Smith PetscCall(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,Mx, My, Mz, 80d0609cedSBarry Smith PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,1, 1,NULL, NULL, NULL,&da3D)); 819566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da3D)); 829566063dSJacob Faibussowitsch PetscCall(DMSetUp(da3D)); 8387adc40fSSajid Ali 8487adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 8587adc40fSSajid Ali Create the parent vector 8687adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/ 879566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da3D, &vec_full)); 889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) vec_full, "full_vector")); 8987adc40fSSajid Ali 9087adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 9187adc40fSSajid Ali Populate the 3D vector 9287adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 939566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da3D, &ixs, &iys, &izs, &ixm, &iym, &izm)); 949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3D, vec_full, &vecdata3d)); 9587adc40fSSajid Ali for (k=izs; k<izs+izm; k++) { 9687adc40fSSajid Ali for (j=iys; j<iys+iym; j++) { 9787adc40fSSajid Ali for (i=ixs; i<ixs+ixm; i++) { 9887adc40fSSajid Ali vecdata3d[k][j][i] = ((i-Mx/2.0)*(j+Mx/2.0))+k*100; 9987adc40fSSajid Ali } 10087adc40fSSajid Ali } 10187adc40fSSajid Ali } 1029566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3D, vec_full, &vecdata3d)); 10387adc40fSSajid Ali 10487adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 10587adc40fSSajid Ali Get an IS corresponding to a 2D slice 10687adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 10787adc40fSSajid Ali if (sliceaxis==DM_X) { 10887adc40fSSajid Ali lower.i = gp; lower.j = 0; lower.k = 0; 10987adc40fSSajid Ali upper.i = gp; upper.j = My; upper.k = Mz; 11087adc40fSSajid Ali } else if (sliceaxis==DM_Y) { 11187adc40fSSajid Ali lower.i = 0; lower.j = gp; lower.k = 0; 11287adc40fSSajid Ali upper.i = Mx; upper.j = gp; upper.k = Mz; 11387adc40fSSajid Ali } else if (sliceaxis==DM_Z) { 11487adc40fSSajid Ali lower.i = 0; lower.j = 0; lower.k = gp; 11587adc40fSSajid Ali upper.i = Mx; upper.j = My; upper.k = gp; 11687adc40fSSajid Ali } 1179566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da3D, &lower, &upper, &patchis_3d, patchis_offproc)); 1189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n IS to select slice from 3D DMDA vector : \n")); 1199566063dSJacob Faibussowitsch PetscCall(ISView(patchis_3d, PETSC_VIEWER_STDOUT_WORLD)); 12087adc40fSSajid Ali 12187adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 12287adc40fSSajid Ali Use the obtained IS to extract the slice as a subvector 12387adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vec_full, patchis_3d, &vec_extracted)); 12587adc40fSSajid Ali 12687adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 12787adc40fSSajid Ali View the extracted subvector 12887adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE)); 1309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in DMDA ordering : \n")); 1319566063dSJacob Faibussowitsch PetscCall(VecView(vec_extracted, PETSC_VIEWER_STDOUT_WORLD)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 13387adc40fSSajid Ali 13487adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 13587adc40fSSajid Ali Query 3D DMDA layout, get the subset MPI communicator 13687adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/ 13787adc40fSSajid Ali if (sliceaxis==DM_X) { 1389566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL)); 1399566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da3D, NULL, &l1, &l2)); 14087adc40fSSajid Ali M1 = My; M2 = Mz; 1419566063dSJacob Faibussowitsch PetscCall(DMDAGetProcessorSubset(da3D, DM_X, gp, &subset_mpi_comm)); 14287adc40fSSajid Ali } else if (sliceaxis==DM_Y) { 1439566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, NULL, &m2, NULL, NULL, NULL, NULL, NULL, NULL)); 1449566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da3D, &l1, NULL, &l2)); 14587adc40fSSajid Ali M1 = Mx; M2 = Mz; 1469566063dSJacob Faibussowitsch PetscCall(DMDAGetProcessorSubset(da3D, DM_Y, gp, &subset_mpi_comm)); 14787adc40fSSajid Ali } else if (sliceaxis==DM_Z) { 1489566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da3D, &l1, &l2, NULL)); 15087adc40fSSajid Ali M1 = Mx; M2 = My; 1519566063dSJacob Faibussowitsch PetscCall(DMDAGetProcessorSubset(da3D, DM_Z, gp, &subset_mpi_comm)); 15287adc40fSSajid Ali } 15387adc40fSSajid Ali 15487adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 15587adc40fSSajid Ali Create 2D DMDA object, 15687adc40fSSajid Ali vector (that will hold the slice as a column major flattened array) & 15787adc40fSSajid Ali index set (that will be used for scattering to the column major 15887adc40fSSajid Ali indexed slice vector) 15987adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/ 16087adc40fSSajid Ali if (subset_mpi_comm != MPI_COMM_NULL) { 1619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subset_mpi_comm, &size)); 1629566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(subset_mpi_comm, "subset MPI subcomm size is : %d, includes global rank : %d \n", size, rank)); 1639566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(subset_mpi_comm, PETSC_STDOUT)); 164d0609cedSBarry Smith PetscCall(DMDACreate2d(subset_mpi_comm,DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M1, M2,m1, m2,1, 1,l1, l2,&da2D)); 1659566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da2D)); 1669566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2D)); 16787adc40fSSajid Ali 16887adc40fSSajid Ali /* Create a 2D patch IS for the slice */ 16987adc40fSSajid Ali lower.i = 0; lower.j = 0; 17087adc40fSSajid Ali upper.i = M1; upper.j = M2; 1719566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da2D, &lower, &upper, &patchis_2d, patchis_offproc)); 17287adc40fSSajid Ali 17387adc40fSSajid Ali /* Convert the 2D patch IS to natural indexing (column major flattened) */ 1749566063dSJacob Faibussowitsch PetscCall(ISDuplicate(patchis_2d, &scatis_natural_slice)); 1759566063dSJacob Faibussowitsch PetscCall(DMDAGetAO(da2D, &da2D_ao)); 1769566063dSJacob Faibussowitsch PetscCall(AOPetscToApplicationIS(da2D_ao, scatis_natural_slice)); 1779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(scatis_natural_slice, &is_array)); 1789566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(scatis_natural_slice, &in)); 17987adc40fSSajid Ali 18087adc40fSSajid Ali /* Create an aliased IS on the 3D DMDA's communicator */ 1819566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, in, is_array, PETSC_USE_POINTER, &scatis_natural_slice_g)); 1829566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(scatis_natural_slice, &is_array)); 18387adc40fSSajid Ali 18487adc40fSSajid Ali /* Create a 2D DMDA global vector */ 1859566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da2D, &vec_slice)); 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) vec_slice, "slice_vector_natural")); 1879566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec_slice ,&vn)); 1889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vec_slice, &array)); 18987adc40fSSajid Ali 19087adc40fSSajid Ali /* Create an aliased version of the above on the 3D DMDA's communicator */ 1919566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, vn, M1*M2, array, &vec_slice_g)); 1929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vec_slice, &array)); 19387adc40fSSajid Ali } else { 19487adc40fSSajid Ali /* Ranks not part of the subset MPI communicator provide no entries, but the routines for creating 19587adc40fSSajid Ali the IS and Vec on the 3D DMDA's communicator still need to called, since they are collective routines */ 1969566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 0, NULL, PETSC_USE_POINTER, &scatis_natural_slice_g)); 1979566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, 0, M1*M2, NULL, &vec_slice_g)); 19887adc40fSSajid Ali } 1999566063dSJacob Faibussowitsch PetscCall(PetscBarrier(NULL)); 20087adc40fSSajid Ali 20187adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 20287adc40fSSajid Ali Create IS that maps from the extracted slice vector 20387adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2049566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vec_extracted, &low, &high)); 2059566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, high-low, low, 1, &scatis_extracted_slice)); 20687adc40fSSajid Ali 20787adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 20887adc40fSSajid Ali Scatter extracted subvector -> natural 2D slice vector 20987adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2109566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec_extracted, scatis_extracted_slice, vec_slice_g, scatis_natural_slice_g, &vscat)); 2119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD)); 2129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD)); 21387adc40fSSajid Ali 21487adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 21587adc40fSSajid Ali View the natural 2D slice vector 21687adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2179566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE)); 2189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in natural ordering : \n")); 2199566063dSJacob Faibussowitsch PetscCall(VecView(vec_slice_g, PETSC_VIEWER_STDOUT_WORLD)); 2209566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 22187adc40fSSajid Ali 22287adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 22387adc40fSSajid Ali Restore subvector 22487adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2259566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vec_full, patchis_3d, &vec_extracted)); 22687adc40fSSajid Ali 22787adc40fSSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 22887adc40fSSajid Ali Destroy data structures and exit. 22987adc40fSSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_full)); 2319566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat)); 2329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&scatis_extracted_slice)); 2339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&scatis_natural_slice_g)); 2349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_slice_g)); 2359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&patchis_3d)); 2369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da3D)); 23787adc40fSSajid Ali 23887adc40fSSajid Ali if (subset_mpi_comm != MPI_COMM_NULL) { 2399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&patchis_2d)); 2409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&scatis_natural_slice)); 2419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec_slice)); 2429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da2D)); 2439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subset_mpi_comm)); 24487adc40fSSajid Ali } 24587adc40fSSajid Ali 2469566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 247b122ec5aSJacob Faibussowitsch return 0; 24887adc40fSSajid Ali } 24987adc40fSSajid Ali 25087adc40fSSajid Ali /*TEST 25187adc40fSSajid Ali 25287adc40fSSajid Ali test: 25387adc40fSSajid Ali nsize: 1 25487adc40fSSajid Ali args: -sliceaxis X -gp 0 25587adc40fSSajid Ali 25687adc40fSSajid Ali test: 25787adc40fSSajid Ali suffix: 2 25887adc40fSSajid Ali nsize: 2 25987adc40fSSajid Ali args: -sliceaxis Y -gp 1 2602071a92dSBarry Smith filter: grep -v "subset MPI subcomm" 26187adc40fSSajid Ali 26287adc40fSSajid Ali test: 26387adc40fSSajid Ali suffix: 3 26487adc40fSSajid Ali nsize: 3 26587adc40fSSajid Ali args: -sliceaxis Z -gp 2 2662071a92dSBarry Smith filter: grep -v "subset MPI subcomm" 26787adc40fSSajid Ali 26887adc40fSSajid Ali test: 26987adc40fSSajid Ali suffix: 4 27087adc40fSSajid Ali nsize: 4 27187adc40fSSajid Ali args: -sliceaxis X -gp 2 2722071a92dSBarry Smith filter: grep -v "subset MPI subcomm" 27387adc40fSSajid Ali 27487adc40fSSajid Ali test: 27587adc40fSSajid Ali suffix: 5 27687adc40fSSajid Ali nsize: 4 27787adc40fSSajid Ali args: -sliceaxis Z -gp 1 2782071a92dSBarry Smith filter: grep -v "subset MPI subcomm" 27987adc40fSSajid Ali 28087adc40fSSajid Ali TEST*/ 281