xref: /petsc/src/dm/tests/ex53.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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\
53b5e53cdSSajid Ali sliceid - set the location where the slice will be extraced 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 
153b5e53cdSSajid Ali int main(int argc,char **argv)
163b5e53cdSSajid Ali {
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   PetscErrorCode ierr;                          /* error checking */
343b5e53cdSSajid Ali 
353b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
363b5e53cdSSajid Ali      Initialize program and set problem parameters
373b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
395f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
405f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
413b5e53cdSSajid Ali 
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, NULL));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-my", &my, NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-sliceid", &sliceid, NULL));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-sliceaxis", &sliceaxis, NULL));
483b5e53cdSSajid Ali 
493b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
503b5e53cdSSajid Ali      Create DMDA object.
513b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
523b5e53cdSSajid Ali   if (dim==3) {
533b5e53cdSSajid Ali     ierr = DMDACreate3d(PETSC_COMM_WORLD,
543b5e53cdSSajid Ali                         DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
553b5e53cdSSajid Ali                         DMDA_STENCIL_STAR,
563b5e53cdSSajid Ali                         mx, my, mz,
573b5e53cdSSajid Ali                         PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,
583b5e53cdSSajid Ali                         1, 1,
593b5e53cdSSajid Ali                         NULL, NULL, NULL,
603b5e53cdSSajid Ali                         &da);CHKERRQ(ierr);
615f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetFromOptions(da));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetUp(da));
633b5e53cdSSajid Ali   } else {
643b5e53cdSSajid Ali     ierr = DMDACreate2d(PETSC_COMM_WORLD,
653b5e53cdSSajid Ali                         DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
663b5e53cdSSajid Ali                         DMDA_STENCIL_STAR,
673b5e53cdSSajid Ali                         mx, my,
683b5e53cdSSajid Ali                         PETSC_DECIDE, PETSC_DECIDE,
693b5e53cdSSajid Ali                         1, 1,
703b5e53cdSSajid Ali                         NULL, NULL,
713b5e53cdSSajid Ali                         &da);CHKERRQ(ierr);
725f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetFromOptions(da));
735f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetUp(da));
743b5e53cdSSajid Ali   }
753b5e53cdSSajid Ali 
763b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
773b5e53cdSSajid Ali      Create the parent vector
783b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da, &vec_full));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) vec_full, "full_vector"));
813b5e53cdSSajid Ali 
823b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
833b5e53cdSSajid Ali      Populate the 3D vector
843b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da, &ixs, &iys, &izs, &ixm, &iym, &izm));
863b5e53cdSSajid Ali   if (dim==3) {
875f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da, vec_full, &vecdata3d));
883b5e53cdSSajid Ali     for (k=izs; k<izs+izm; k++) {
893b5e53cdSSajid Ali       for (j=iys; j<iys+iym; j++) {
903b5e53cdSSajid Ali         for (i=ixs; i<ixs+ixm; i++) {
913b5e53cdSSajid Ali           vecdata3d[k][j][i] = ((i-mx/2)*(j+mx/2))+k*100;
923b5e53cdSSajid Ali         }
933b5e53cdSSajid Ali       }
943b5e53cdSSajid Ali     }
955f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da, vec_full, &vecdata3d));
963b5e53cdSSajid Ali   } else {
975f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da, vec_full, &vecdata2d));
983b5e53cdSSajid Ali     for (j=iys; j<iys+iym; j++) {
993b5e53cdSSajid Ali       for (i=ixs; i<ixs+ixm; i++) {
1003b5e53cdSSajid Ali         vecdata2d[j][i] = ((i-mx/2)*(j+mx/2));
1013b5e53cdSSajid Ali       }
1023b5e53cdSSajid Ali     }
1035f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da, vec_full, &vecdata2d));
1043b5e53cdSSajid Ali   }
1053b5e53cdSSajid Ali 
1063b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1073b5e53cdSSajid Ali      Get an IS corresponding to a 2D slice
1083b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1093b5e53cdSSajid Ali   if (sliceaxis==0) {
1103b5e53cdSSajid Ali     lower.i = sliceid; lower.j = 0;  lower.k = 0;  lower.c = 1;
1113b5e53cdSSajid Ali     upper.i = sliceid; upper.j = my; upper.k = mz; upper.c = 1;
1123b5e53cdSSajid Ali   } else if (sliceaxis==1) {
1133b5e53cdSSajid Ali     lower.i = 0;  lower.j = sliceid; lower.k = 0;  lower.c = 1;
1143b5e53cdSSajid Ali     upper.i = mx; upper.j = sliceid; upper.k = mz; upper.c = 1;
1153b5e53cdSSajid Ali   } else if (sliceaxis==2 && dim==3) {
1163b5e53cdSSajid Ali     lower.i = 0;  lower.j = 0;  lower.k = sliceid; lower.c = 1;
1173b5e53cdSSajid Ali     upper.i = mx; upper.j = my; upper.k = sliceid; upper.c = 1;
1183b5e53cdSSajid Ali   }
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreatePatchIS(da, &lower, &upper, &selectis, patchis_offproc));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(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      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSubVector(vec_full, selectis, &vec_slice));
1263b5e53cdSSajid Ali 
1273b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1283b5e53cdSSajid Ali      View the extracted subvector
1293b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(vec_slice, PETSC_VIEWER_STDOUT_WORLD));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
1333b5e53cdSSajid Ali 
1343b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1353b5e53cdSSajid Ali      Restore subvector, destroy data structures and exit.
1363b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreSubVector(vec_full, selectis, &vec_slice));
1383b5e53cdSSajid Ali 
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&selectis));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vec_full));
1423b5e53cdSSajid Ali 
143*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
144*b122ec5aSJacob 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