xref: /petsc/src/dm/tests/ex53.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
343b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
353b5e53cdSSajid Ali      Initialize program and set problem parameters
363b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
403b5e53cdSSajid Ali 
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, NULL));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, NULL));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-sliceid", &sliceid, NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-sliceaxis", &sliceaxis, NULL));
473b5e53cdSSajid Ali 
483b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
493b5e53cdSSajid Ali      Create DMDA object.
503b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
513b5e53cdSSajid Ali   if (dim==3) {
52*d0609cedSBarry Smith     PetscCall(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,mx, my, mz,
53*d0609cedSBarry Smith                            PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,1, 1,NULL, NULL, NULL,&da));
549566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(da));
559566063dSJacob Faibussowitsch     PetscCall(DMSetUp(da));
563b5e53cdSSajid Ali   } else {
57*d0609cedSBarry Smith     PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,mx, my,PETSC_DECIDE, PETSC_DECIDE,
58*d0609cedSBarry Smith                            1, 1,NULL, NULL,&da));
599566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(da));
609566063dSJacob Faibussowitsch     PetscCall(DMSetUp(da));
613b5e53cdSSajid Ali   }
623b5e53cdSSajid Ali 
633b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
643b5e53cdSSajid Ali      Create the parent vector
653b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
669566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &vec_full));
679566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) vec_full, "full_vector"));
683b5e53cdSSajid Ali 
693b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
703b5e53cdSSajid Ali      Populate the 3D vector
713b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
729566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &ixs, &iys, &izs, &ixm, &iym, &izm));
733b5e53cdSSajid Ali   if (dim==3) {
749566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, vec_full, &vecdata3d));
753b5e53cdSSajid Ali     for (k=izs; k<izs+izm; k++) {
763b5e53cdSSajid Ali       for (j=iys; j<iys+iym; j++) {
773b5e53cdSSajid Ali         for (i=ixs; i<ixs+ixm; i++) {
783b5e53cdSSajid Ali           vecdata3d[k][j][i] = ((i-mx/2)*(j+mx/2))+k*100;
793b5e53cdSSajid Ali         }
803b5e53cdSSajid Ali       }
813b5e53cdSSajid Ali     }
829566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, vec_full, &vecdata3d));
833b5e53cdSSajid Ali   } else {
849566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, vec_full, &vecdata2d));
853b5e53cdSSajid Ali     for (j=iys; j<iys+iym; j++) {
863b5e53cdSSajid Ali       for (i=ixs; i<ixs+ixm; i++) {
873b5e53cdSSajid Ali         vecdata2d[j][i] = ((i-mx/2)*(j+mx/2));
883b5e53cdSSajid Ali       }
893b5e53cdSSajid Ali     }
909566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, vec_full, &vecdata2d));
913b5e53cdSSajid Ali   }
923b5e53cdSSajid Ali 
933b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
943b5e53cdSSajid Ali      Get an IS corresponding to a 2D slice
953b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
963b5e53cdSSajid Ali   if (sliceaxis==0) {
973b5e53cdSSajid Ali     lower.i = sliceid; lower.j = 0;  lower.k = 0;  lower.c = 1;
983b5e53cdSSajid Ali     upper.i = sliceid; upper.j = my; upper.k = mz; upper.c = 1;
993b5e53cdSSajid Ali   } else if (sliceaxis==1) {
1003b5e53cdSSajid Ali     lower.i = 0;  lower.j = sliceid; lower.k = 0;  lower.c = 1;
1013b5e53cdSSajid Ali     upper.i = mx; upper.j = sliceid; upper.k = mz; upper.c = 1;
1023b5e53cdSSajid Ali   } else if (sliceaxis==2 && dim==3) {
1033b5e53cdSSajid Ali     lower.i = 0;  lower.j = 0;  lower.k = sliceid; lower.c = 1;
1043b5e53cdSSajid Ali     upper.i = mx; upper.j = my; upper.k = sliceid; upper.c = 1;
1053b5e53cdSSajid Ali   }
1069566063dSJacob Faibussowitsch   PetscCall(DMDACreatePatchIS(da, &lower, &upper, &selectis, patchis_offproc));
1079566063dSJacob Faibussowitsch   PetscCall(ISView(selectis, PETSC_VIEWER_STDOUT_WORLD));
1083b5e53cdSSajid Ali 
1093b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1103b5e53cdSSajid Ali      Use the obtained IS to extract the slice as a subvector
1113b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1129566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(vec_full, selectis, &vec_slice));
1133b5e53cdSSajid Ali 
1143b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1153b5e53cdSSajid Ali      View the extracted subvector
1163b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1179566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE));
1189566063dSJacob Faibussowitsch   PetscCall(VecView(vec_slice, PETSC_VIEWER_STDOUT_WORLD));
1199566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
1203b5e53cdSSajid Ali 
1213b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1223b5e53cdSSajid Ali      Restore subvector, destroy data structures and exit.
1233b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1249566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(vec_full, selectis, &vec_slice));
1253b5e53cdSSajid Ali 
1269566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selectis));
1279566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vec_full));
1293b5e53cdSSajid Ali 
1309566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
131b122ec5aSJacob Faibussowitsch   return 0;
1323b5e53cdSSajid Ali }
1333b5e53cdSSajid Ali 
1343b5e53cdSSajid Ali /*TEST
1353b5e53cdSSajid Ali 
1363b5e53cdSSajid Ali     test:
1373b5e53cdSSajid Ali       nsize: 1
1383b5e53cdSSajid Ali       args: -sliceaxis 0
1393b5e53cdSSajid Ali 
1403b5e53cdSSajid Ali     test:
1413b5e53cdSSajid Ali       suffix: 2
1423b5e53cdSSajid Ali       nsize:  2
1433b5e53cdSSajid Ali       args: -sliceaxis 1
1443b5e53cdSSajid Ali 
1453b5e53cdSSajid Ali     test:
1463b5e53cdSSajid Ali       suffix: 3
1473b5e53cdSSajid Ali       nsize:  4
1483b5e53cdSSajid Ali       args:  -sliceaxis 2
1493b5e53cdSSajid Ali 
1503b5e53cdSSajid Ali     test:
1513b5e53cdSSajid Ali       suffix: 4
1523b5e53cdSSajid Ali       nsize:  2
1533b5e53cdSSajid Ali       args: -sliceaxis 1 -dim 2
1543b5e53cdSSajid Ali 
1553b5e53cdSSajid Ali     test:
1563b5e53cdSSajid Ali       suffix: 5
1573b5e53cdSSajid Ali       nsize:  3
1583b5e53cdSSajid Ali       args: -sliceaxis 0 -dim 2
1593b5e53cdSSajid Ali 
1603b5e53cdSSajid Ali TEST*/
161