xref: /petsc/src/dm/tutorials/ex22.c (revision 87adc40fb667840c167c5c478289a4942dce160c)
1*87adc40fSSajid Ali static char help[] ="Extract a 2D slice in natural ordering from a 3D vector, Command line options : \n\
2*87adc40fSSajid Ali Mx/My/Mz - set the dimensions of the parent vector \n\
3*87adc40fSSajid Ali sliceaxis - string describing the axis along which the sice will be selected : X, Y, Z \n\
4*87adc40fSSajid Ali gp - global grid point number along the sliceaxis direction where the slice will be extracted from the parent vector \n";
5*87adc40fSSajid Ali 
6*87adc40fSSajid Ali /*
7*87adc40fSSajid Ali   This example shows to extract a 2D slice in natural ordering
8*87adc40fSSajid Ali   from a 3D DMDA vector (first by extracting the slice and then
9*87adc40fSSajid Ali   by converting it to natural ordering)
10*87adc40fSSajid Ali */
11*87adc40fSSajid Ali 
12*87adc40fSSajid Ali #include <petscdmda.h>
13*87adc40fSSajid Ali 
14*87adc40fSSajid Ali const char *const sliceaxes[] = {"X","Y","Z","sliceaxis","DM_",NULL};
15*87adc40fSSajid Ali 
16*87adc40fSSajid Ali int main(int argc,char **argv)
17*87adc40fSSajid Ali {
18*87adc40fSSajid Ali   DM                da3D;                          /* 3D DMDA object */
19*87adc40fSSajid Ali   DM                da2D;                          /* 2D DMDA object */
20*87adc40fSSajid Ali   Vec               vec_full;                      /* Parent vector */
21*87adc40fSSajid Ali   Vec               vec_extracted;                 /* Extracted slice vector (in DMDA ordering) */
22*87adc40fSSajid Ali   Vec               vec_slice;                     /* vector in natural ordering */
23*87adc40fSSajid Ali   Vec               vec_slice_g;                   /* aliased vector in natural ordering */
24*87adc40fSSajid Ali   IS                patchis_3d;                    /* IS to select slice and extract subvector */
25*87adc40fSSajid Ali   IS                patchis_2d;                    /* Patch IS for 2D vector, will be converted to application ordering */
26*87adc40fSSajid Ali   IS                scatis_extracted_slice;        /* PETSc indexed IS for extracted slice */
27*87adc40fSSajid Ali   IS                scatis_natural_slice;          /* natural/application ordered IS for slice*/
28*87adc40fSSajid Ali   IS                scatis_natural_slice_g;        /* aliased natural/application ordered IS  for slice */
29*87adc40fSSajid Ali   VecScatter        vscat;                         /* scatter slice in DMDA ordering <-> slice in column major ordering */
30*87adc40fSSajid Ali   AO                da2D_ao;                       /* AO associated with 2D DMDA */
31*87adc40fSSajid Ali   MPI_Comm          subset_mpi_comm=MPI_COMM_NULL; /* MPI communicator where the slice lives */
32*87adc40fSSajid Ali   PetscScalar       ***vecdata3d;                  /* Pointer to access 3d parent vector */
33*87adc40fSSajid Ali   const PetscScalar *array;                        /* pointer to create aliased Vec */
34*87adc40fSSajid Ali   PetscInt          Mx=4,My=4,Mz=4;                /* Dimensions for 3D DMDA */
35*87adc40fSSajid Ali   const PetscInt    *l1,*l2;                       /* 3D DMDA layout */
36*87adc40fSSajid Ali   PetscInt          M1=-1,M2=-1;                   /* Dimensions for 2D DMDA */
37*87adc40fSSajid Ali   PetscInt          m1=-1,m2=-1;                   /* Layouts for 2D DMDA */
38*87adc40fSSajid Ali   PetscInt          gp=2;                          /* grid point along sliceaxis to pick the slice */
39*87adc40fSSajid Ali   DMDirection       sliceaxis=DM_X;                /* Select axis along which the slice will be extracted */
40*87adc40fSSajid Ali   PetscInt          i,j,k;                         /* Iteration indices */
41*87adc40fSSajid Ali   PetscInt          ixs,iys,izs;                   /* Corner indices for 3D vector */
42*87adc40fSSajid Ali   PetscInt          ixm,iym,izm;                   /* Widths of parent vector */
43*87adc40fSSajid Ali   PetscInt          low, high;                     /* ownership range indices */
44*87adc40fSSajid Ali   PetscInt          in;                            /* local size index for IS*/
45*87adc40fSSajid Ali   PetscInt          vn;                            /* local size index */
46*87adc40fSSajid Ali   const PetscInt    *is_array;                     /* pointer to create aliased IS */
47*87adc40fSSajid Ali   MatStencil        lower, upper;                  /* Stencils to select slice for Vec */
48*87adc40fSSajid Ali   PetscBool         patchis_offproc = PETSC_FALSE; /* flag to DMDACreatePatchIS indicating that off-proc values are to be ignored */
49*87adc40fSSajid Ali   PetscMPIInt       rank,size;                     /* MPI rank and size */
50*87adc40fSSajid Ali   PetscErrorCode    ierr;                          /* error checking */
51*87adc40fSSajid Ali 
52*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53*87adc40fSSajid Ali      Initialize program and set problem parameters
54*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55*87adc40fSSajid Ali   ierr = PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr;
56*87adc40fSSajid Ali   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
57*87adc40fSSajid Ali   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
58*87adc40fSSajid Ali 
59*87adc40fSSajid Ali   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "ex22 DMDA tutorial example options", "DMDA");CHKERRQ(ierr);
60*87adc40fSSajid Ali   ierr = PetscOptionsRangeInt("-Mx", "dimension along x-axis", "ex22.c", Mx, &Mx, NULL, 0, PETSC_MAX_INT);CHKERRQ(ierr);
61*87adc40fSSajid Ali   ierr = PetscOptionsRangeInt("-My", "dimension along y-axis", "ex22.c", My, &My, NULL, 0, PETSC_MAX_INT);CHKERRQ(ierr);
62*87adc40fSSajid Ali   ierr = PetscOptionsRangeInt("-Mz", "dimension along z-axis", "ex22.c", Mz, &Mz, NULL, 0, PETSC_MAX_INT);CHKERRQ(ierr);
63*87adc40fSSajid Ali   ierr = PetscOptionsEnum("-sliceaxis","axis along which 2D slice is extracted from : X, Y, Z","",sliceaxes,(PetscEnum)sliceaxis,(PetscEnum*)&sliceaxis,NULL);CHKERRQ(ierr);
64*87adc40fSSajid Ali   ierr = PetscOptionsRangeInt("-gp", "index along sliceaxis at which 2D slice is extracted", "ex22.c", gp, &gp, NULL, 0, PETSC_MAX_INT);CHKERRQ(ierr);
65*87adc40fSSajid Ali   ierr = PetscOptionsEnd();CHKERRQ(ierr);
66*87adc40fSSajid Ali 
67*87adc40fSSajid Ali   /* Ensure that the requested slice is not out of bounds for the selected axis */
68*87adc40fSSajid Ali   if (sliceaxis==DM_X) {
69*87adc40fSSajid Ali     if (gp>Mx) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "grid point along sliceaxis is larger than largest index!");
70*87adc40fSSajid Ali   } else if (sliceaxis==DM_Y) {
71*87adc40fSSajid Ali     if (gp>My) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "grid point along sliceaxis is larger than largest index!");
72*87adc40fSSajid Ali   } else if (sliceaxis==DM_Z) {
73*87adc40fSSajid Ali     if (gp>Mz) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "grid point along sliceaxis is larger than largest index!");
74*87adc40fSSajid Ali   }
75*87adc40fSSajid Ali 
76*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77*87adc40fSSajid Ali      Create 3D DMDA object.
78*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
79*87adc40fSSajid Ali   ierr = DMDACreate3d(PETSC_COMM_WORLD,
80*87adc40fSSajid Ali                       DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
81*87adc40fSSajid Ali                       DMDA_STENCIL_STAR,
82*87adc40fSSajid Ali                       Mx, My, Mz,
83*87adc40fSSajid Ali                       PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,
84*87adc40fSSajid Ali                       1, 1,
85*87adc40fSSajid Ali                       NULL, NULL, NULL,
86*87adc40fSSajid Ali                       &da3D);CHKERRQ(ierr);
87*87adc40fSSajid Ali   ierr = DMSetFromOptions(da3D);CHKERRQ(ierr);
88*87adc40fSSajid Ali   ierr = DMSetUp(da3D);CHKERRQ(ierr);
89*87adc40fSSajid Ali 
90*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91*87adc40fSSajid Ali      Create the parent vector
92*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
93*87adc40fSSajid Ali   ierr = DMCreateGlobalVector(da3D, &vec_full);CHKERRQ(ierr);
94*87adc40fSSajid Ali   ierr = PetscObjectSetName((PetscObject) vec_full, "full_vector");CHKERRQ(ierr);
95*87adc40fSSajid Ali 
96*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97*87adc40fSSajid Ali      Populate the 3D vector
98*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99*87adc40fSSajid Ali   ierr = DMDAGetCorners(da3D, &ixs, &iys, &izs, &ixm, &iym, &izm);CHKERRQ(ierr);
100*87adc40fSSajid Ali   ierr = DMDAVecGetArray(da3D, vec_full, &vecdata3d);CHKERRQ(ierr);
101*87adc40fSSajid Ali   for (k=izs; k<izs+izm; k++){
102*87adc40fSSajid Ali     for (j=iys; j<iys+iym; j++) {
103*87adc40fSSajid Ali       for (i=ixs; i<ixs+ixm; i++) {
104*87adc40fSSajid Ali         vecdata3d[k][j][i] = ((i-Mx/2.0)*(j+Mx/2.0))+k*100;
105*87adc40fSSajid Ali       }
106*87adc40fSSajid Ali     }
107*87adc40fSSajid Ali   }
108*87adc40fSSajid Ali   ierr = DMDAVecRestoreArray(da3D, vec_full, &vecdata3d);CHKERRQ(ierr);
109*87adc40fSSajid Ali 
110*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111*87adc40fSSajid Ali      Get an IS corresponding to a 2D slice
112*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113*87adc40fSSajid Ali   if (sliceaxis==DM_X) {
114*87adc40fSSajid Ali     lower.i = gp; lower.j = 0;  lower.k = 0;
115*87adc40fSSajid Ali     upper.i = gp; upper.j = My; upper.k = Mz;
116*87adc40fSSajid Ali   } else if (sliceaxis==DM_Y) {
117*87adc40fSSajid Ali     lower.i = 0;  lower.j = gp; lower.k = 0;
118*87adc40fSSajid Ali     upper.i = Mx; upper.j = gp; upper.k = Mz;
119*87adc40fSSajid Ali   } else if (sliceaxis==DM_Z) {
120*87adc40fSSajid Ali     lower.i = 0;  lower.j = 0;  lower.k = gp;
121*87adc40fSSajid Ali     upper.i = Mx; upper.j = My; upper.k = gp;
122*87adc40fSSajid Ali   }
123*87adc40fSSajid Ali   ierr = DMDACreatePatchIS(da3D, &lower, &upper, &patchis_3d, patchis_offproc);CHKERRQ(ierr);
124*87adc40fSSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD, "\n IS to select slice from 3D DMDA vector : \n");CHKERRQ(ierr);
125*87adc40fSSajid Ali   ierr = ISView(patchis_3d, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
126*87adc40fSSajid Ali 
127*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128*87adc40fSSajid Ali      Use the obtained IS to extract the slice as a subvector
129*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130*87adc40fSSajid Ali   ierr = VecGetSubVector(vec_full, patchis_3d, &vec_extracted);CHKERRQ(ierr);
131*87adc40fSSajid Ali 
132*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133*87adc40fSSajid Ali      View the extracted subvector
134*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135*87adc40fSSajid Ali   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE);CHKERRQ(ierr);
136*87adc40fSSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in DMDA ordering : \n");CHKERRQ(ierr);
137*87adc40fSSajid Ali   ierr = VecView(vec_extracted, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
138*87adc40fSSajid Ali   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
139*87adc40fSSajid Ali 
140*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141*87adc40fSSajid Ali      Query 3D DMDA layout, get the subset MPI communicator
142*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
143*87adc40fSSajid Ali   if (sliceaxis==DM_X) {
144*87adc40fSSajid Ali     ierr = DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
145*87adc40fSSajid Ali     ierr = DMDAGetOwnershipRanges(da3D, NULL, &l1, &l2);CHKERRQ(ierr);
146*87adc40fSSajid Ali     M1 = My; M2 = Mz;
147*87adc40fSSajid Ali     ierr = DMDAGetProcessorSubset(da3D, DM_X, gp, &subset_mpi_comm);CHKERRQ(ierr);
148*87adc40fSSajid Ali   } else if (sliceaxis==DM_Y) {
149*87adc40fSSajid Ali     ierr = DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, NULL, &m2, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
150*87adc40fSSajid Ali     ierr = DMDAGetOwnershipRanges(da3D, &l1, NULL, &l2);CHKERRQ(ierr);
151*87adc40fSSajid Ali     M1 = Mx; M2 = Mz;
152*87adc40fSSajid Ali     ierr = DMDAGetProcessorSubset(da3D, DM_Y, gp, &subset_mpi_comm);CHKERRQ(ierr);
153*87adc40fSSajid Ali   } else if (sliceaxis==DM_Z) {
154*87adc40fSSajid Ali     ierr = DMDAGetInfo(da3D, NULL, NULL, NULL, NULL, &m1, &m2, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
155*87adc40fSSajid Ali     ierr = DMDAGetOwnershipRanges(da3D, &l1, &l2, NULL);CHKERRQ(ierr);
156*87adc40fSSajid Ali     M1 = Mx; M2 = My;
157*87adc40fSSajid Ali     ierr = DMDAGetProcessorSubset(da3D, DM_Z, gp, &subset_mpi_comm);CHKERRQ(ierr);
158*87adc40fSSajid Ali   }
159*87adc40fSSajid Ali 
160*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161*87adc40fSSajid Ali      Create 2D DMDA object,
162*87adc40fSSajid Ali      vector (that will hold the slice as a column major flattened array) &
163*87adc40fSSajid Ali      index set (that will be used for scattering to the column major
164*87adc40fSSajid Ali      indexed slice vector)
165*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
166*87adc40fSSajid Ali   if (subset_mpi_comm != MPI_COMM_NULL) {
167*87adc40fSSajid Ali     ierr = MPI_Comm_size(subset_mpi_comm, &size);CHKERRMPI(ierr);
168*87adc40fSSajid Ali     ierr = PetscSynchronizedPrintf(subset_mpi_comm, "subset MPI subcomm size is : %d, includes global rank : %d \n", size, rank);CHKERRQ(ierr);
169*87adc40fSSajid Ali     ierr = PetscSynchronizedFlush(subset_mpi_comm, PETSC_STDOUT);CHKERRQ(ierr);
170*87adc40fSSajid Ali     ierr = DMDACreate2d(subset_mpi_comm,
171*87adc40fSSajid Ali                         DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
172*87adc40fSSajid Ali                         DMDA_STENCIL_STAR,
173*87adc40fSSajid Ali                         M1, M2,
174*87adc40fSSajid Ali                         m1, m2,
175*87adc40fSSajid Ali                         1, 1,
176*87adc40fSSajid Ali                         l1, l2,
177*87adc40fSSajid Ali                         &da2D);CHKERRQ(ierr);
178*87adc40fSSajid Ali     ierr = DMSetFromOptions(da2D);CHKERRQ(ierr);
179*87adc40fSSajid Ali     ierr = DMSetUp(da2D);CHKERRQ(ierr);
180*87adc40fSSajid Ali 
181*87adc40fSSajid Ali     /* Create a 2D patch IS for the slice */
182*87adc40fSSajid Ali     lower.i = 0;  lower.j = 0;
183*87adc40fSSajid Ali     upper.i = M1; upper.j = M2;
184*87adc40fSSajid Ali     ierr = DMDACreatePatchIS(da2D, &lower, &upper, &patchis_2d, patchis_offproc);CHKERRQ(ierr);
185*87adc40fSSajid Ali 
186*87adc40fSSajid Ali     /* Convert the 2D patch IS to natural indexing (column major flattened) */
187*87adc40fSSajid Ali     ierr = ISDuplicate(patchis_2d, &scatis_natural_slice);CHKERRQ(ierr);
188*87adc40fSSajid Ali     ierr = DMDAGetAO(da2D, &da2D_ao);CHKERRQ(ierr);
189*87adc40fSSajid Ali     ierr = AOPetscToApplicationIS(da2D_ao, scatis_natural_slice);CHKERRQ(ierr);
190*87adc40fSSajid Ali     ierr = ISGetIndices(scatis_natural_slice, &is_array);CHKERRQ(ierr);
191*87adc40fSSajid Ali     ierr = ISGetLocalSize(scatis_natural_slice, &in);CHKERRQ(ierr);
192*87adc40fSSajid Ali 
193*87adc40fSSajid Ali     /* Create an aliased IS on the 3D DMDA's communicator */
194*87adc40fSSajid Ali     ierr = ISCreateGeneral(PETSC_COMM_WORLD, in, is_array, PETSC_USE_POINTER, &scatis_natural_slice_g);CHKERRQ(ierr);
195*87adc40fSSajid Ali     ierr = ISRestoreIndices(scatis_natural_slice, &is_array);CHKERRQ(ierr);
196*87adc40fSSajid Ali 
197*87adc40fSSajid Ali     /* Create a 2D DMDA global vector */
198*87adc40fSSajid Ali     ierr = DMCreateGlobalVector(da2D, &vec_slice);CHKERRQ(ierr);
199*87adc40fSSajid Ali     ierr = PetscObjectSetName((PetscObject) vec_slice, "slice_vector_natural");CHKERRQ(ierr);
200*87adc40fSSajid Ali     ierr = VecGetLocalSize(vec_slice ,&vn);CHKERRQ(ierr);
201*87adc40fSSajid Ali     ierr = VecGetArrayRead(vec_slice, &array);CHKERRQ(ierr);
202*87adc40fSSajid Ali 
203*87adc40fSSajid Ali     /* Create an aliased version of the above on the 3D DMDA's communicator */
204*87adc40fSSajid Ali     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, vn, M1*M2, array, &vec_slice_g);CHKERRQ(ierr);
205*87adc40fSSajid Ali     ierr = VecRestoreArrayRead(vec_slice, &array);CHKERRQ(ierr);
206*87adc40fSSajid Ali   } else {
207*87adc40fSSajid Ali     /* Ranks not part of the subset MPI communicator provide no entries, but the routines for creating
208*87adc40fSSajid Ali        the IS and Vec on the 3D DMDA's communicator still need to called, since they are collective routines */
209*87adc40fSSajid Ali     ierr = ISCreateGeneral(PETSC_COMM_WORLD, 0, NULL, PETSC_USE_POINTER, &scatis_natural_slice_g);CHKERRQ(ierr);
210*87adc40fSSajid Ali     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, 0, M1*M2, NULL, &vec_slice_g);CHKERRQ(ierr);
211*87adc40fSSajid Ali   }
212*87adc40fSSajid Ali   ierr = PetscBarrier(NULL);CHKERRQ(ierr);
213*87adc40fSSajid Ali 
214*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215*87adc40fSSajid Ali      Create IS that maps from the extracted slice vector
216*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217*87adc40fSSajid Ali   ierr = VecGetOwnershipRange(vec_extracted, &low, &high);CHKERRQ(ierr);
218*87adc40fSSajid Ali   ierr = ISCreateStride(PETSC_COMM_WORLD, high-low, low, 1, &scatis_extracted_slice);CHKERRQ(ierr);
219*87adc40fSSajid Ali 
220*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221*87adc40fSSajid Ali      Scatter extracted subvector -> natural 2D slice vector
222*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223*87adc40fSSajid Ali   ierr = VecScatterCreate(vec_extracted, scatis_extracted_slice, vec_slice_g, scatis_natural_slice_g, &vscat);CHKERRQ(ierr);
224*87adc40fSSajid Ali   ierr = VecScatterBegin(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
225*87adc40fSSajid Ali   ierr = VecScatterEnd(vscat, vec_extracted, vec_slice_g, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
226*87adc40fSSajid Ali 
227*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228*87adc40fSSajid Ali      View the natural 2D slice vector
229*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230*87adc40fSSajid Ali   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE);CHKERRQ(ierr);
231*87adc40fSSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD, "\n Extracted slice vector, in natural ordering : \n");CHKERRQ(ierr);
232*87adc40fSSajid Ali   ierr = VecView(vec_slice_g, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
233*87adc40fSSajid Ali   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
234*87adc40fSSajid Ali 
235*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236*87adc40fSSajid Ali      Restore subvector
237*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238*87adc40fSSajid Ali   ierr = VecRestoreSubVector(vec_full, patchis_3d, &vec_extracted);CHKERRQ(ierr);
239*87adc40fSSajid Ali 
240*87adc40fSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241*87adc40fSSajid Ali      Destroy data structures and exit.
242*87adc40fSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243*87adc40fSSajid Ali   ierr = VecDestroy(&vec_full);CHKERRQ(ierr);
244*87adc40fSSajid Ali   ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
245*87adc40fSSajid Ali   ierr = ISDestroy(&scatis_extracted_slice);CHKERRQ(ierr);
246*87adc40fSSajid Ali   ierr = ISDestroy(&scatis_natural_slice_g);CHKERRQ(ierr);
247*87adc40fSSajid Ali   ierr = VecDestroy(&vec_slice_g);CHKERRQ(ierr);
248*87adc40fSSajid Ali   ierr = ISDestroy(&patchis_3d);CHKERRQ(ierr);
249*87adc40fSSajid Ali   ierr = DMDestroy(&da3D);CHKERRQ(ierr);
250*87adc40fSSajid Ali 
251*87adc40fSSajid Ali   if (subset_mpi_comm != MPI_COMM_NULL) {
252*87adc40fSSajid Ali     ierr = ISDestroy(&patchis_2d);CHKERRQ(ierr);
253*87adc40fSSajid Ali     ierr = ISDestroy(&scatis_natural_slice);CHKERRQ(ierr);
254*87adc40fSSajid Ali     ierr = VecDestroy(&vec_slice);CHKERRQ(ierr);
255*87adc40fSSajid Ali     ierr = DMDestroy(&da2D);CHKERRQ(ierr);
256*87adc40fSSajid Ali     ierr = MPI_Comm_free(&subset_mpi_comm);CHKERRMPI(ierr);
257*87adc40fSSajid Ali   }
258*87adc40fSSajid Ali 
259*87adc40fSSajid Ali   ierr = PetscFinalize();
260*87adc40fSSajid Ali   return ierr;
261*87adc40fSSajid Ali }
262*87adc40fSSajid Ali 
263*87adc40fSSajid Ali /*TEST
264*87adc40fSSajid Ali 
265*87adc40fSSajid Ali     test:
266*87adc40fSSajid Ali       nsize: 1
267*87adc40fSSajid Ali       args: -sliceaxis X -gp 0
268*87adc40fSSajid Ali 
269*87adc40fSSajid Ali     test:
270*87adc40fSSajid Ali       suffix: 2
271*87adc40fSSajid Ali       nsize:  2
272*87adc40fSSajid Ali       args: -sliceaxis Y -gp 1
273*87adc40fSSajid Ali 
274*87adc40fSSajid Ali     test:
275*87adc40fSSajid Ali       suffix: 3
276*87adc40fSSajid Ali       nsize:  3
277*87adc40fSSajid Ali       args:  -sliceaxis Z -gp 2
278*87adc40fSSajid Ali 
279*87adc40fSSajid Ali     test:
280*87adc40fSSajid Ali       suffix: 4
281*87adc40fSSajid Ali       nsize:  4
282*87adc40fSSajid Ali       args: -sliceaxis X -gp 2
283*87adc40fSSajid Ali 
284*87adc40fSSajid Ali     test:
285*87adc40fSSajid Ali       suffix: 5
286*87adc40fSSajid Ali       nsize:  4
287*87adc40fSSajid Ali       args: -sliceaxis Z -gp 1
288*87adc40fSSajid Ali 
289*87adc40fSSajid Ali TEST*/
290