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