xref: /petsc/src/mat/tests/ex211.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() in parallel.";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown PetscErrorCode ISGetSeqIS_SameColDist_Private(Mat mat,IS isrow,IS iscol,IS *isrow_d,IS *iscol_d,IS *iscol_o,const PetscInt *garray[])
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscErrorCode ierr;
10*c4762a1bSJed Brown   Vec            x,cmap;
11*c4762a1bSJed Brown   const PetscInt *is_idx;
12*c4762a1bSJed Brown   PetscScalar    *xarray,*cmaparray;
13*c4762a1bSJed Brown   PetscInt       ncols,isstart,*idx,m,rstart,count;
14*c4762a1bSJed Brown   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)mat->data;
15*c4762a1bSJed Brown   Mat            B=a->B;
16*c4762a1bSJed Brown   Vec            lvec=a->lvec,lcmap;
17*c4762a1bSJed Brown   PetscInt       i,cstart,cend,Bn=B->cmap->N;
18*c4762a1bSJed Brown   MPI_Comm       comm;
19*c4762a1bSJed Brown   PetscMPIInt    rank;
20*c4762a1bSJed Brown   VecScatter     Mvctx;
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   PetscFunctionBegin;
23*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   /* (1) iscol is a sub-column vector of mat, pad it with '-1.' to form a full vector x */
28*c4762a1bSJed Brown   ierr = MatCreateVecs(mat,&x,NULL);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = VecDuplicate(x,&cmap);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = VecSet(x,-1.0);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = VecSet(cmap,-1.0);CHKERRQ(ierr);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   ierr = VecDuplicate(lvec,&lcmap);CHKERRQ(ierr);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   /* Get start indices */
36*c4762a1bSJed Brown   ierr = MPI_Scan(&ncols,&isstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
37*c4762a1bSJed Brown   isstart -= ncols;
38*c4762a1bSJed Brown   ierr = MatGetOwnershipRangeColumn(mat,&cstart,&cend);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   ierr = ISGetIndices(iscol,&is_idx);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = VecGetArray(cmap,&cmaparray);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscMalloc1(ncols,&idx);CHKERRQ(ierr);
44*c4762a1bSJed Brown   for (i=0; i<ncols; i++) {
45*c4762a1bSJed Brown     xarray[is_idx[i]-cstart]    = (PetscScalar)is_idx[i];
46*c4762a1bSJed Brown     cmaparray[is_idx[i]-cstart] = (PetscScalar)(i + isstart);      /* global index of iscol[i] */
47*c4762a1bSJed Brown     idx[i]                      = is_idx[i]-cstart; /* local index of iscol[i]  */
48*c4762a1bSJed Brown   }
49*c4762a1bSJed Brown   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = VecRestoreArray(cmap,&cmaparray);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscol,&is_idx);CHKERRQ(ierr);
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown   /* Get iscol_d */
54*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,iscol_d);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = ISGetBlockSize(iscol,&i);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = ISSetBlockSize(*iscol_d,i);CHKERRQ(ierr);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   /* Get isrow_d */
59*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
60*c4762a1bSJed Brown   rstart = mat->rmap->rstart;
61*c4762a1bSJed Brown   ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = ISGetIndices(isrow,&is_idx);CHKERRQ(ierr);
63*c4762a1bSJed Brown   for (i=0; i<m; i++) idx[i] = is_idx[i]-rstart;
64*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrow,&is_idx);CHKERRQ(ierr);
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,isrow_d);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = ISGetBlockSize(isrow,&i);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = ISSetBlockSize(*isrow_d,i);CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   /* (2) Scatter x and cmap using aij->Mvctx to get their off-process portions (see MatMult_MPIAIJ) */
71*c4762a1bSJed Brown #if 0
72*c4762a1bSJed Brown   if (!a->Mvctx_mpi1) {
73*c4762a1bSJed Brown     /* a->Mvctx causes random 'count' in o-build? See src/mat/tests/runex59_2 */
74*c4762a1bSJed Brown     a->Mvctx_mpi1_flg = PETSC_TRUE;
75*c4762a1bSJed Brown     ierr = MatSetUpMultiply_MPIAIJ(mat);CHKERRQ(ierr);
76*c4762a1bSJed Brown   }
77*c4762a1bSJed Brown   Mvctx = a->Mvctx_mpi1;
78*c4762a1bSJed Brown #endif
79*c4762a1bSJed Brown   Mvctx = a->Mvctx;
80*c4762a1bSJed Brown   ierr = VecScatterBegin(Mvctx,x,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = VecScatterEnd(Mvctx,x,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown   ierr = VecScatterBegin(Mvctx,cmap,lcmap,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = VecScatterEnd(Mvctx,cmap,lcmap,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   /* (3) create sequential iscol_o (a subset of iscol) and isgarray */
87*c4762a1bSJed Brown   /* off-process column indices */
88*c4762a1bSJed Brown   count = 0;
89*c4762a1bSJed Brown   PetscInt *cmap1;
90*c4762a1bSJed Brown   ierr = PetscMalloc1(Bn,&idx);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = PetscMalloc1(Bn,&cmap1);CHKERRQ(ierr);
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   ierr = VecGetArray(lvec,&xarray);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = VecGetArray(lcmap,&cmaparray);CHKERRQ(ierr);
95*c4762a1bSJed Brown   for (i=0; i<Bn; i++) {
96*c4762a1bSJed Brown     if (PetscRealPart(xarray[i]) > -1.0) {
97*c4762a1bSJed Brown       idx[count]   = i;                   /* local column index in off-diagonal part B */
98*c4762a1bSJed Brown       cmap1[count] = (PetscInt)(PetscRealPart(cmaparray[i]));  /* column index in submat */
99*c4762a1bSJed Brown       count++;
100*c4762a1bSJed Brown     }
101*c4762a1bSJed Brown   }
102*c4762a1bSJed Brown   printf("[%d] Bn %d, count %d\n",rank,Bn,count);
103*c4762a1bSJed Brown   ierr = VecRestoreArray(lvec,&xarray);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = VecRestoreArray(lcmap,&cmaparray);CHKERRQ(ierr);
105*c4762a1bSJed Brown   if (count != 6) {
106*c4762a1bSJed Brown     printf("[%d] count %d != 6 lvec:\n",rank,count);
107*c4762a1bSJed Brown     ierr = VecView(lvec,0);CHKERRQ(ierr);
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown     printf("[%d] count %d != 6 lcmap:\n",rank,count);
110*c4762a1bSJed Brown     ierr = VecView(lcmap,0);CHKERRQ(ierr);
111*c4762a1bSJed Brown     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"count %d != 6",count);
112*c4762a1bSJed Brown   }
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,idx,PETSC_COPY_VALUES,iscol_o);CHKERRQ(ierr);
115*c4762a1bSJed Brown   /* cannot ensure iscol_o has same blocksize as iscol! */
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown   ierr = PetscFree(idx);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   *garray = cmap1;
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = VecDestroy(&cmap);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = VecDestroy(&lcmap);CHKERRQ(ierr);
124*c4762a1bSJed Brown   PetscFunctionReturn(0);
125*c4762a1bSJed Brown }
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown int main(int argc,char **args)
128*c4762a1bSJed Brown {
129*c4762a1bSJed Brown   Mat            C,A;
130*c4762a1bSJed Brown   PetscInt       i,j,m = 3,n = 2,rstart,rend;
131*c4762a1bSJed Brown   PetscMPIInt    size,rank;
132*c4762a1bSJed Brown   PetscErrorCode ierr;
133*c4762a1bSJed Brown   PetscScalar    v;
134*c4762a1bSJed Brown   IS             isrow,iscol;
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
137*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
139*c4762a1bSJed Brown   n    = 2*size;
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown   /*
147*c4762a1bSJed Brown         This is JUST to generate a nice test matrix, all processors fill up
148*c4762a1bSJed Brown     the entire matrix. This is not something one would ever do in practice.
149*c4762a1bSJed Brown   */
150*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
151*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
152*c4762a1bSJed Brown     for (j=0; j<m*n; j++) {
153*c4762a1bSJed Brown       v    = i + j + 1;
154*c4762a1bSJed Brown       ierr = MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
155*c4762a1bSJed Brown     }
156*c4762a1bSJed Brown   }
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown   /*
162*c4762a1bSJed Brown      Generate a new matrix consisting of every second row and column of
163*c4762a1bSJed Brown    the original matrix
164*c4762a1bSJed Brown   */
165*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
166*c4762a1bSJed Brown   /* Create parallel IS with the rows we want on THIS processor */
167*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&isrow);CHKERRQ(ierr);
168*c4762a1bSJed Brown   /* Create parallel IS with the rows we want on THIS processor (same as rows for now) */
169*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&iscol);CHKERRQ(ierr);
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown   IS             iscol_d,isrow_d,iscol_o;
172*c4762a1bSJed Brown   const PetscInt *garray;
173*c4762a1bSJed Brown   ierr = ISGetSeqIS_SameColDist_Private(C,isrow,iscol,&isrow_d,&iscol_d,&iscol_o,&garray);CHKERRQ(ierr);
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown   ierr = ISDestroy(&isrow_d);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = ISDestroy(&iscol_d);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = ISDestroy(&iscol_o);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = PetscFree(garray);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   ierr = MatCreateSubMatrix(C,isrow,iscol,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = MatCreateSubMatrix(C,isrow,iscol,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown   ierr = ISDestroy(&isrow);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = PetscFinalize();
188*c4762a1bSJed Brown   return ierr;
189*c4762a1bSJed Brown }
190