xref: /petsc/src/dm/tests/ex53.c (revision 3b5e53cdb169e94712b800a2d1c9e10766bebf6e)
1*3b5e53cdSSajid Ali static char help[] ="Use DMDACreatePatchIS  to extract a slice from a vector, Command line options :\n\
2*3b5e53cdSSajid Ali mx/my/mz - set the dimensions of the parent vector\n\
3*3b5e53cdSSajid Ali dim - set the dimensionality of the parent vector (2,3)\n\
4*3b5e53cdSSajid Ali sliceaxis - Integer describing the axis along which the sice will be selected (0-X, 1-Y, 2-Z)\n\
5*3b5e53cdSSajid Ali sliceid - set the location where the slice will be extraced from the parent vector\n";
6*3b5e53cdSSajid Ali 
7*3b5e53cdSSajid Ali /*
8*3b5e53cdSSajid Ali    This test checks the functionality of DMDACreatePatchIS when
9*3b5e53cdSSajid Ali    extracting a 2D vector from a 3D vector and 1D vector from a
10*3b5e53cdSSajid Ali    2D vector.
11*3b5e53cdSSajid Ali    */
12*3b5e53cdSSajid Ali 
13*3b5e53cdSSajid Ali #include <petscdmda.h>
14*3b5e53cdSSajid Ali 
15*3b5e53cdSSajid Ali int main(int argc,char **argv)
16*3b5e53cdSSajid Ali {
17*3b5e53cdSSajid Ali   PetscMPIInt    rank, size;                    /* MPI rank and size */
18*3b5e53cdSSajid Ali   PetscInt       mx=4,my=4,mz=4;                /* Dimensions of parent vector */
19*3b5e53cdSSajid Ali   PetscInt       sliceid=2;                     /* k (z) index to pick the slice */
20*3b5e53cdSSajid Ali   PetscInt       sliceaxis=2;                   /* Select axis along which the slice will be extracted */
21*3b5e53cdSSajid Ali   PetscInt       dim=3;                         /* Dimension of the parent vector */
22*3b5e53cdSSajid Ali   PetscInt       i,j,k;                         /* Iteration indices */
23*3b5e53cdSSajid Ali   PetscInt       ixs,iys,izs;                   /* Corner indices for 3D vector */
24*3b5e53cdSSajid Ali   PetscInt       ixm,iym,izm;                   /* Widths of parent vector */
25*3b5e53cdSSajid Ali   PetscScalar    ***vecdata3d;                  /* Pointer to access 3d parent vector */
26*3b5e53cdSSajid Ali   PetscScalar    **vecdata2d;                   /* Pointer to access 2d parent vector */
27*3b5e53cdSSajid Ali   DM             da;                            /* 2D/3D DMDA object */
28*3b5e53cdSSajid Ali   Vec            vec_full;                      /* Parent vector */
29*3b5e53cdSSajid Ali   Vec            vec_slice;                     /* Slice vector */
30*3b5e53cdSSajid Ali   MatStencil     lower, upper;                  /* Stencils to select slice */
31*3b5e53cdSSajid Ali   IS             selectis;                      /* IS to select slice and extract subvector */
32*3b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_FALSE; /* flag to DMDACreatePatchIS indicating that off-proc values are to be ignored */
33*3b5e53cdSSajid Ali   PetscErrorCode ierr;                          /* error checking */
34*3b5e53cdSSajid Ali 
35*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36*3b5e53cdSSajid Ali      Initialize program and set problem parameters
37*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38*3b5e53cdSSajid Ali   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
39*3b5e53cdSSajid Ali   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
40*3b5e53cdSSajid Ali   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
41*3b5e53cdSSajid Ali 
42*3b5e53cdSSajid Ali   ierr = PetscOptionsGetInt(NULL, NULL, "-mx", &mx, NULL);CHKERRQ(ierr);
43*3b5e53cdSSajid Ali   ierr = PetscOptionsGetInt(NULL, NULL, "-my", &my, NULL);CHKERRQ(ierr);
44*3b5e53cdSSajid Ali   ierr = PetscOptionsGetInt(NULL, NULL, "-mz", &mz, NULL);CHKERRQ(ierr);
45*3b5e53cdSSajid Ali   ierr = PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL);CHKERRQ(ierr);
46*3b5e53cdSSajid Ali   ierr = PetscOptionsGetInt(NULL, NULL, "-sliceid", &sliceid, NULL);CHKERRQ(ierr);
47*3b5e53cdSSajid Ali   ierr = PetscOptionsGetInt(NULL, NULL, "-sliceaxis", &sliceaxis, NULL);CHKERRQ(ierr);
48*3b5e53cdSSajid Ali 
49*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50*3b5e53cdSSajid Ali      Create DMDA object.
51*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
52*3b5e53cdSSajid Ali   if (dim==3) {
53*3b5e53cdSSajid Ali     ierr = DMDACreate3d(PETSC_COMM_WORLD,
54*3b5e53cdSSajid Ali                         DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
55*3b5e53cdSSajid Ali                         DMDA_STENCIL_STAR,
56*3b5e53cdSSajid Ali                         mx, my, mz,
57*3b5e53cdSSajid Ali                         PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,
58*3b5e53cdSSajid Ali                         1, 1,
59*3b5e53cdSSajid Ali                         NULL, NULL, NULL,
60*3b5e53cdSSajid Ali                         &da);CHKERRQ(ierr);
61*3b5e53cdSSajid Ali     ierr = DMSetFromOptions(da);CHKERRQ(ierr);
62*3b5e53cdSSajid Ali     ierr = DMSetUp(da);CHKERRQ(ierr);
63*3b5e53cdSSajid Ali   } else {
64*3b5e53cdSSajid Ali     ierr = DMDACreate2d(PETSC_COMM_WORLD,
65*3b5e53cdSSajid Ali                         DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
66*3b5e53cdSSajid Ali                         DMDA_STENCIL_STAR,
67*3b5e53cdSSajid Ali                         mx, my,
68*3b5e53cdSSajid Ali                         PETSC_DECIDE, PETSC_DECIDE,
69*3b5e53cdSSajid Ali                         1, 1,
70*3b5e53cdSSajid Ali                         NULL, NULL,
71*3b5e53cdSSajid Ali                         &da);CHKERRQ(ierr);
72*3b5e53cdSSajid Ali     ierr = DMSetFromOptions(da);CHKERRQ(ierr);
73*3b5e53cdSSajid Ali     ierr = DMSetUp(da);CHKERRQ(ierr);
74*3b5e53cdSSajid Ali   }
75*3b5e53cdSSajid Ali 
76*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77*3b5e53cdSSajid Ali      Create the parent vector
78*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -*/
79*3b5e53cdSSajid Ali   ierr = DMCreateGlobalVector(da, &vec_full);CHKERRQ(ierr);
80*3b5e53cdSSajid Ali   ierr = PetscObjectSetName((PetscObject) vec_full, "full_vector");CHKERRQ(ierr);
81*3b5e53cdSSajid Ali 
82*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83*3b5e53cdSSajid Ali      Populate the 3D vector
84*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85*3b5e53cdSSajid Ali   ierr = DMDAGetCorners(da, &ixs, &iys, &izs, &ixm, &iym, &izm);CHKERRQ(ierr);
86*3b5e53cdSSajid Ali   if (dim==3) {
87*3b5e53cdSSajid Ali     ierr = DMDAVecGetArray(da, vec_full, &vecdata3d);CHKERRQ(ierr);
88*3b5e53cdSSajid Ali     for (k=izs; k<izs+izm; k++){
89*3b5e53cdSSajid Ali       for (j=iys; j<iys+iym; j++) {
90*3b5e53cdSSajid Ali         for (i=ixs; i<ixs+ixm; i++) {
91*3b5e53cdSSajid Ali           vecdata3d[k][j][i] = ((i-mx/2)*(j+mx/2))+k*100;
92*3b5e53cdSSajid Ali         }
93*3b5e53cdSSajid Ali       }
94*3b5e53cdSSajid Ali     }
95*3b5e53cdSSajid Ali     ierr = DMDAVecRestoreArray(da, vec_full, &vecdata3d);CHKERRQ(ierr);
96*3b5e53cdSSajid Ali   } else {
97*3b5e53cdSSajid Ali     ierr = DMDAVecGetArray(da, vec_full, &vecdata2d);CHKERRQ(ierr);
98*3b5e53cdSSajid Ali     for (j=iys; j<iys+iym; j++) {
99*3b5e53cdSSajid Ali       for (i=ixs; i<ixs+ixm; i++) {
100*3b5e53cdSSajid Ali         vecdata2d[j][i] = ((i-mx/2)*(j+mx/2));
101*3b5e53cdSSajid Ali       }
102*3b5e53cdSSajid Ali     }
103*3b5e53cdSSajid Ali     ierr = DMDAVecRestoreArray(da, vec_full, &vecdata2d);CHKERRQ(ierr);
104*3b5e53cdSSajid Ali   }
105*3b5e53cdSSajid Ali 
106*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107*3b5e53cdSSajid Ali      Get an IS corresponding to a 2D slice
108*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109*3b5e53cdSSajid Ali   if (sliceaxis==0) {
110*3b5e53cdSSajid Ali     lower.i = sliceid; lower.j = 0;  lower.k = 0;  lower.c = 1;
111*3b5e53cdSSajid Ali     upper.i = sliceid; upper.j = my; upper.k = mz; upper.c = 1;
112*3b5e53cdSSajid Ali   } else if (sliceaxis==1) {
113*3b5e53cdSSajid Ali     lower.i = 0;  lower.j = sliceid; lower.k = 0;  lower.c = 1;
114*3b5e53cdSSajid Ali     upper.i = mx; upper.j = sliceid; upper.k = mz; upper.c = 1;
115*3b5e53cdSSajid Ali   } else if (sliceaxis==2 && dim==3) {
116*3b5e53cdSSajid Ali     lower.i = 0;  lower.j = 0;  lower.k = sliceid; lower.c = 1;
117*3b5e53cdSSajid Ali     upper.i = mx; upper.j = my; upper.k = sliceid; upper.c = 1;
118*3b5e53cdSSajid Ali   }
119*3b5e53cdSSajid Ali   ierr = DMDACreatePatchIS(da, &lower, &upper, &selectis, patchis_offproc);CHKERRQ(ierr);
120*3b5e53cdSSajid Ali   ierr = ISView(selectis, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
121*3b5e53cdSSajid Ali 
122*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123*3b5e53cdSSajid Ali      Use the obtained IS to extract the slice as a subvector
124*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125*3b5e53cdSSajid Ali   ierr = VecGetSubVector(vec_full, selectis, &vec_slice);CHKERRQ(ierr);
126*3b5e53cdSSajid Ali 
127*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128*3b5e53cdSSajid Ali      View the extracted subvector
129*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130*3b5e53cdSSajid Ali   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE);CHKERRQ(ierr);
131*3b5e53cdSSajid Ali   ierr = VecView(vec_slice, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
132*3b5e53cdSSajid Ali   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
133*3b5e53cdSSajid Ali 
134*3b5e53cdSSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135*3b5e53cdSSajid Ali      Restore subvector, destroy data structures and exit.
136*3b5e53cdSSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137*3b5e53cdSSajid Ali   ierr = VecRestoreSubVector(vec_full, selectis, &vec_slice);CHKERRQ(ierr);
138*3b5e53cdSSajid Ali 
139*3b5e53cdSSajid Ali   ierr = ISDestroy(&selectis);CHKERRQ(ierr);
140*3b5e53cdSSajid Ali   ierr = DMDestroy(&da);CHKERRQ(ierr);
141*3b5e53cdSSajid Ali   ierr = VecDestroy(&vec_full);CHKERRQ(ierr);
142*3b5e53cdSSajid Ali 
143*3b5e53cdSSajid Ali   ierr = PetscFinalize();
144*3b5e53cdSSajid Ali   return ierr;
145*3b5e53cdSSajid Ali }
146*3b5e53cdSSajid Ali 
147*3b5e53cdSSajid Ali /*TEST
148*3b5e53cdSSajid Ali 
149*3b5e53cdSSajid Ali     test:
150*3b5e53cdSSajid Ali       nsize: 1
151*3b5e53cdSSajid Ali       args: -sliceaxis 0
152*3b5e53cdSSajid Ali 
153*3b5e53cdSSajid Ali     test:
154*3b5e53cdSSajid Ali       suffix: 2
155*3b5e53cdSSajid Ali       nsize:  2
156*3b5e53cdSSajid Ali       args: -sliceaxis 1
157*3b5e53cdSSajid Ali 
158*3b5e53cdSSajid Ali     test:
159*3b5e53cdSSajid Ali       suffix: 3
160*3b5e53cdSSajid Ali       nsize:  4
161*3b5e53cdSSajid Ali       args:  -sliceaxis 2
162*3b5e53cdSSajid Ali 
163*3b5e53cdSSajid Ali     test:
164*3b5e53cdSSajid Ali       suffix: 4
165*3b5e53cdSSajid Ali       nsize:  2
166*3b5e53cdSSajid Ali       args: -sliceaxis 1 -dim 2
167*3b5e53cdSSajid Ali 
168*3b5e53cdSSajid Ali     test:
169*3b5e53cdSSajid Ali       suffix: 5
170*3b5e53cdSSajid Ali       nsize:  3
171*3b5e53cdSSajid Ali       args: -sliceaxis 0 -dim 2
172*3b5e53cdSSajid Ali 
173*3b5e53cdSSajid Ali TEST*/
174