xref: /petsc/src/dm/tutorials/ex22.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode    ierr;                          /* error checking */
5187adc40fSSajid Ali 
5287adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5387adc40fSSajid Ali      Initialize program and set problem parameters
5487adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, (char*)0, help));
565f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
575f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
5887adc40fSSajid Ali 
5987adc40fSSajid Ali   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "ex22 DMDA tutorial example options", "DMDA");CHKERRQ(ierr);
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsRangeInt("-Mx", "dimension along x-axis", "ex22.c", Mx, &Mx, NULL, 0, PETSC_MAX_INT));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsRangeInt("-My", "dimension along y-axis", "ex22.c", My, &My, NULL, 0, PETSC_MAX_INT));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsRangeInt("-Mz", "dimension along z-axis", "ex22.c", Mz, &Mz, NULL, 0, PETSC_MAX_INT));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-sliceaxis","axis along which 2D slice is extracted from : X, Y, Z","",sliceaxes,(PetscEnum)sliceaxis,(PetscEnum*)&sliceaxis,NULL));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsRangeInt("-gp", "index along sliceaxis at which 2D slice is extracted", "ex22.c", gp, &gp, NULL, 0, PETSC_MAX_INT));
6587adc40fSSajid Ali   ierr = PetscOptionsEnd();CHKERRQ(ierr);
6687adc40fSSajid Ali 
6787adc40fSSajid Ali   /* Ensure that the requested slice is not out of bounds for the selected axis */
6887adc40fSSajid Ali   if (sliceaxis==DM_X) {
692c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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) {
712c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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) {
732c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
7987adc40fSSajid Ali   ierr = DMDACreate3d(PETSC_COMM_WORLD,
8087adc40fSSajid Ali                       DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
8187adc40fSSajid Ali                       DMDA_STENCIL_STAR,
8287adc40fSSajid Ali                       Mx, My, Mz,
8387adc40fSSajid Ali                       PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,
8487adc40fSSajid Ali                       1, 1,
8587adc40fSSajid Ali                       NULL, NULL, NULL,
8687adc40fSSajid Ali                       &da3D);CHKERRQ(ierr);
875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da3D));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da3D));
8987adc40fSSajid Ali 
9087adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9187adc40fSSajid Ali      Create the parent vector
9287adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da3D, &vec_full));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) vec_full, "full_vector"));
9587adc40fSSajid Ali 
9687adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9787adc40fSSajid Ali      Populate the 3D vector
9887adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da3D, &ixs, &iys, &izs, &ixm, &iym, &izm));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da3D, vec_full, &vecdata3d));
10187adc40fSSajid Ali   for (k=izs; k<izs+izm; k++) {
10287adc40fSSajid Ali     for (j=iys; j<iys+iym; j++) {
10387adc40fSSajid Ali       for (i=ixs; i<ixs+ixm; i++) {
10487adc40fSSajid Ali         vecdata3d[k][j][i] = ((i-Mx/2.0)*(j+Mx/2.0))+k*100;
10587adc40fSSajid Ali       }
10687adc40fSSajid Ali     }
10787adc40fSSajid Ali   }
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da3D, vec_full, &vecdata3d));
10987adc40fSSajid Ali 
11087adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11187adc40fSSajid Ali      Get an IS corresponding to a 2D slice
11287adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
11387adc40fSSajid Ali   if (sliceaxis==DM_X) {
11487adc40fSSajid Ali     lower.i = gp; lower.j = 0;  lower.k = 0;
11587adc40fSSajid Ali     upper.i = gp; upper.j = My; upper.k = Mz;
11687adc40fSSajid Ali   } else if (sliceaxis==DM_Y) {
11787adc40fSSajid Ali     lower.i = 0;  lower.j = gp; lower.k = 0;
11887adc40fSSajid Ali     upper.i = Mx; upper.j = gp; upper.k = Mz;
11987adc40fSSajid Ali   } else if (sliceaxis==DM_Z) {
12087adc40fSSajid Ali     lower.i = 0;  lower.j = 0;  lower.k = gp;
12187adc40fSSajid Ali     upper.i = Mx; upper.j = My; upper.k = gp;
12287adc40fSSajid Ali   }
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreatePatchIS(da3D, &lower, &upper, &patchis_3d, patchis_offproc));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\n IS to select slice from 3D DMDA vector : \n"));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(patchis_3d, PETSC_VIEWER_STDOUT_WORLD));
12687adc40fSSajid Ali 
12787adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
12887adc40fSSajid Ali      Use the obtained IS to extract the slice as a subvector
12987adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSubVector(vec_full, patchis_3d, &vec_extracted));
13187adc40fSSajid Ali 
13287adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
13387adc40fSSajid Ali      View the extracted subvector
13487adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in DMDA ordering : \n"));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(vec_extracted, PETSC_VIEWER_STDOUT_WORLD));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
13987adc40fSSajid Ali 
14087adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14187adc40fSSajid Ali      Query 3D DMDA layout, get the subset MPI communicator
14287adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
14387adc40fSSajid Ali   if (sliceaxis==DM_X) {
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL));
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetOwnershipRanges(da3D, NULL, &l1, &l2));
14687adc40fSSajid Ali     M1 = My; M2 = Mz;
1475f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetProcessorSubset(da3D, DM_X, gp, &subset_mpi_comm));
14887adc40fSSajid Ali   } else if (sliceaxis==DM_Y) {
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, NULL, &m2, NULL, NULL, NULL, NULL, NULL, NULL));
1505f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetOwnershipRanges(da3D, &l1, NULL, &l2));
15187adc40fSSajid Ali     M1 = Mx; M2 = Mz;
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetProcessorSubset(da3D, DM_Y, gp, &subset_mpi_comm));
15387adc40fSSajid Ali   } else if (sliceaxis==DM_Z) {
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetOwnershipRanges(da3D, &l1, &l2, NULL));
15687adc40fSSajid Ali     M1 = Mx; M2 = My;
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetProcessorSubset(da3D, DM_Z, gp, &subset_mpi_comm));
15887adc40fSSajid Ali   }
15987adc40fSSajid Ali 
16087adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
16187adc40fSSajid Ali      Create 2D DMDA object,
16287adc40fSSajid Ali      vector (that will hold the slice as a column major flattened array) &
16387adc40fSSajid Ali      index set (that will be used for scattering to the column major
16487adc40fSSajid Ali      indexed slice vector)
16587adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
16687adc40fSSajid Ali   if (subset_mpi_comm != MPI_COMM_NULL) {
1675f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_size(subset_mpi_comm, &size));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedPrintf(subset_mpi_comm, "subset MPI subcomm size is : %d, includes global rank : %d \n", size, rank));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFlush(subset_mpi_comm, PETSC_STDOUT));
17087adc40fSSajid Ali     ierr = DMDACreate2d(subset_mpi_comm,
17187adc40fSSajid Ali                         DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
17287adc40fSSajid Ali                         DMDA_STENCIL_STAR,
17387adc40fSSajid Ali                         M1, M2,
17487adc40fSSajid Ali                         m1, m2,
17587adc40fSSajid Ali                         1, 1,
17687adc40fSSajid Ali                         l1, l2,
17787adc40fSSajid Ali                         &da2D);CHKERRQ(ierr);
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetFromOptions(da2D));
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetUp(da2D));
18087adc40fSSajid Ali 
18187adc40fSSajid Ali     /* Create a 2D patch IS for the slice */
18287adc40fSSajid Ali     lower.i = 0;  lower.j = 0;
18387adc40fSSajid Ali     upper.i = M1; upper.j = M2;
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreatePatchIS(da2D, &lower, &upper, &patchis_2d, patchis_offproc));
18587adc40fSSajid Ali 
18687adc40fSSajid Ali     /* Convert the 2D patch IS to natural indexing (column major flattened) */
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDuplicate(patchis_2d, &scatis_natural_slice));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetAO(da2D, &da2D_ao));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(AOPetscToApplicationIS(da2D_ao, scatis_natural_slice));
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(scatis_natural_slice, &is_array));
1915f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(scatis_natural_slice, &in));
19287adc40fSSajid Ali 
19387adc40fSSajid Ali     /* Create an aliased IS on the 3D DMDA's communicator */
1945f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD, in, is_array, PETSC_USE_POINTER, &scatis_natural_slice_g));
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(scatis_natural_slice, &is_array));
19687adc40fSSajid Ali 
19787adc40fSSajid Ali     /* Create a 2D DMDA global vector */
1985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateGlobalVector(da2D, &vec_slice));
1995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) vec_slice, "slice_vector_natural"));
2005f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(vec_slice ,&vn));
2015f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(vec_slice, &array));
20287adc40fSSajid Ali 
20387adc40fSSajid Ali     /* Create an aliased version of the above on the 3D DMDA's communicator */
2045f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, vn, M1*M2, array, &vec_slice_g));
2055f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(vec_slice, &array));
20687adc40fSSajid Ali   } else {
20787adc40fSSajid Ali     /* Ranks not part of the subset MPI communicator provide no entries, but the routines for creating
20887adc40fSSajid Ali        the IS and Vec on the 3D DMDA's communicator still need to called, since they are collective routines */
2095f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD, 0, NULL, PETSC_USE_POINTER, &scatis_natural_slice_g));
2105f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, 0, M1*M2, NULL, &vec_slice_g));
21187adc40fSSajid Ali   }
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBarrier(NULL));
21387adc40fSSajid Ali 
21487adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
21587adc40fSSajid Ali      Create IS that maps from the extracted slice vector
21687adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(vec_extracted, &low, &high));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD, high-low, low, 1, &scatis_extracted_slice));
21987adc40fSSajid Ali 
22087adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22187adc40fSSajid Ali      Scatter extracted subvector -> natural 2D slice vector
22287adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(vec_extracted, scatis_extracted_slice, vec_slice_g, scatis_natural_slice_g, &vscat));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD));
22687adc40fSSajid Ali 
22787adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22887adc40fSSajid Ali      View the natural 2D slice vector
22987adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE));
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in natural ordering : \n"));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(vec_slice_g, PETSC_VIEWER_STDOUT_WORLD));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
23487adc40fSSajid Ali 
23587adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23687adc40fSSajid Ali      Restore subvector
23787adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreSubVector(vec_full, patchis_3d, &vec_extracted));
23987adc40fSSajid Ali 
24087adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24187adc40fSSajid Ali      Destroy data structures and exit.
24287adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vec_full));
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&scatis_extracted_slice));
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&scatis_natural_slice_g));
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vec_slice_g));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&patchis_3d));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da3D));
25087adc40fSSajid Ali 
25187adc40fSSajid Ali   if (subset_mpi_comm != MPI_COMM_NULL) {
2525f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&patchis_2d));
2535f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&scatis_natural_slice));
2545f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&vec_slice));
2555f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&da2D));
2565f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_free(&subset_mpi_comm));
25787adc40fSSajid Ali   }
25887adc40fSSajid Ali 
259*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
260*b122ec5aSJacob Faibussowitsch   return 0;
26187adc40fSSajid Ali }
26287adc40fSSajid Ali 
26387adc40fSSajid Ali /*TEST
26487adc40fSSajid Ali 
26587adc40fSSajid Ali     test:
26687adc40fSSajid Ali       nsize: 1
26787adc40fSSajid Ali       args: -sliceaxis X -gp 0
26887adc40fSSajid Ali 
26987adc40fSSajid Ali     test:
27087adc40fSSajid Ali       suffix: 2
27187adc40fSSajid Ali       nsize:  2
27287adc40fSSajid Ali       args: -sliceaxis Y -gp 1
2732071a92dSBarry Smith       filter: grep -v "subset MPI subcomm"
27487adc40fSSajid Ali 
27587adc40fSSajid Ali     test:
27687adc40fSSajid Ali       suffix: 3
27787adc40fSSajid Ali       nsize:  3
27887adc40fSSajid Ali       args:  -sliceaxis Z -gp 2
2792071a92dSBarry Smith       filter: grep -v "subset MPI subcomm"
28087adc40fSSajid Ali 
28187adc40fSSajid Ali     test:
28287adc40fSSajid Ali       suffix: 4
28387adc40fSSajid Ali       nsize:  4
28487adc40fSSajid Ali       args: -sliceaxis X -gp 2
2852071a92dSBarry Smith       filter: grep -v "subset MPI subcomm"
28687adc40fSSajid Ali 
28787adc40fSSajid Ali     test:
28887adc40fSSajid Ali       suffix: 5
28987adc40fSSajid Ali       nsize:  4
29087adc40fSSajid Ali       args: -sliceaxis Z -gp 1
2912071a92dSBarry Smith       filter: grep -v "subset MPI subcomm"
29287adc40fSSajid Ali 
29387adc40fSSajid Ali TEST*/
294