xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 1b53664090dae02bf12c6ba2f5a1edc60516922c)
1be1d678aSKris Buschelman 
2b97cf49bSBarry Smith /*
3b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
4b97cf49bSBarry Smith */
5e37d9405SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h>    /*I "petscmat.h" I*/
6841d17a1SFande Kong #include <petscsf.h>
7b97cf49bSBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
10dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
11b97cf49bSBarry Smith {
123eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
13dfbe8321SBarry Smith   PetscErrorCode    ierr;
14d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n;
152dcb1b2aSMatthew Knepley   const char        *name;
16ce1411ecSBarry Smith   PetscViewerFormat format;
17b97cf49bSBarry Smith 
18433994e6SBarry Smith   PetscFunctionBegin;
193a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
20b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
21fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
223a40ed3dSBarry Smith     PetscFunctionReturn(0);
23fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
24ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported");
250752156aSBarry Smith   } else {
26d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
277b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
28b97cf49bSBarry Smith     for (i=0; i<m; i++) {
29d0f46423SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr);
30b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
3177431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
320752156aSBarry Smith       }
33b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
34b97cf49bSBarry Smith     }
35d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
36b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
377b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
387b23a99aSBarry Smith   }
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
40b97cf49bSBarry Smith }
41b97cf49bSBarry Smith 
424a2ae208SSatish Balay #undef __FUNCT__
434a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
44dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
45b97cf49bSBarry Smith {
46dfbe8321SBarry Smith   PetscErrorCode ierr;
47ace3abfcSBarry Smith   PetscBool      iascii;
48b97cf49bSBarry Smith 
49433994e6SBarry Smith   PetscFunctionBegin;
50251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5132077d6dSBarry Smith   if (iascii) {
523eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
53b97cf49bSBarry Smith   }
543a40ed3dSBarry Smith   PetscFunctionReturn(0);
55b97cf49bSBarry Smith }
56b97cf49bSBarry Smith 
574a2ae208SSatish Balay #undef __FUNCT__
584a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
59dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat)
60b97cf49bSBarry Smith {
613eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
62dfbe8321SBarry Smith   PetscErrorCode ierr;
6394d884c6SBarry Smith 
6494d884c6SBarry Smith   PetscFunctionBegin;
65aa482453SBarry Smith #if defined(PETSC_USE_LOG)
66d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
67b97cf49bSBarry Smith #endif
6805b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
690400557aSBarry Smith   if (a->freeaij) {
7014196391SBarry Smith     if (a->freeaijwithfree) {
7114196391SBarry Smith       if (a->i) free(a->i);
7214196391SBarry Smith       if (a->j) free(a->j);
7314196391SBarry Smith     } else {
74606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
75606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
7605b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
7714196391SBarry Smith     }
780400557aSBarry Smith   }
792b1d8763SJed Brown   ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
80bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
81dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
82bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
83bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
843a40ed3dSBarry Smith   PetscFunctionReturn(0);
85b97cf49bSBarry Smith }
86b97cf49bSBarry Smith 
874a2ae208SSatish Balay #undef __FUNCT__
884a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
89ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
90b97cf49bSBarry Smith {
913eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
9263ba0a88SBarry Smith   PetscErrorCode ierr;
93b97cf49bSBarry Smith 
94433994e6SBarry Smith   PetscFunctionBegin;
9512c028f9SKris Buschelman   switch (op) {
9677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9712c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
989a4540c5SBarry Smith   case MAT_HERMITIAN:
994e0d8c25SBarry Smith     a->symmetric = flg;
1009a4540c5SBarry Smith     break;
1019a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1029a4540c5SBarry Smith     break;
10312c028f9SKris Buschelman   default:
104290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
10512c028f9SKris Buschelman     break;
106b97cf49bSBarry Smith   }
1073a40ed3dSBarry Smith   PetscFunctionReturn(0);
108b97cf49bSBarry Smith }
109b97cf49bSBarry Smith 
110b97cf49bSBarry Smith 
111b97cf49bSBarry Smith /*
112b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
113b97cf49bSBarry Smith */
114b97cf49bSBarry Smith 
1154a2ae208SSatish Balay #undef __FUNCT__
1164a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
117dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
118b97cf49bSBarry Smith {
1193eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1206849ba73SBarry Smith   PetscErrorCode ierr;
121d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n;
122b97cf49bSBarry Smith 
123433994e6SBarry Smith   PetscFunctionBegin;
124785e854fSJed Brown   ierr = PetscMalloc1(m,&a->diag);CHKERRQ(ierr);
1253bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,m*sizeof(PetscInt));CHKERRQ(ierr);
126d0f46423SBarry Smith   for (i=0; i<A->rmap->n; i++) {
127b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
128b97cf49bSBarry Smith       if (a->j[j] == i) {
12909f38230SBarry Smith         a->diag[i] = j;
130b97cf49bSBarry Smith         break;
131b97cf49bSBarry Smith       }
132b97cf49bSBarry Smith     }
133b97cf49bSBarry Smith   }
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
135b97cf49bSBarry Smith }
136b97cf49bSBarry Smith 
1374a2ae208SSatish Balay #undef __FUNCT__
1384a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
139b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
140b97cf49bSBarry Smith {
1413eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
1422b1d8763SJed Brown   PetscErrorCode ierr;
143b97cf49bSBarry Smith 
144433994e6SBarry Smith   PetscFunctionBegin;
145d0f46423SBarry Smith   row -= A->rmap->rstart;
1460752156aSBarry Smith 
147e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
148b97cf49bSBarry Smith 
149b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
1502b1d8763SJed Brown   if (v) {
1512b1d8763SJed Brown     PetscInt j;
1522b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
1532b1d8763SJed Brown       ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
1542b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz);
155785e854fSJed Brown       ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr);
156b97cf49bSBarry Smith     }
1572b1d8763SJed Brown     for (j=0; j<*nz; j++) a->rowvalues[j] = a->values[a->i[row]+j];
1582b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
159b97cf49bSBarry Smith   }
16092b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
1613a40ed3dSBarry Smith   PetscFunctionReturn(0);
162b97cf49bSBarry Smith }
163b97cf49bSBarry Smith 
1644a2ae208SSatish Balay #undef __FUNCT__
1654a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
166b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
167b97cf49bSBarry Smith {
168433994e6SBarry Smith   PetscFunctionBegin;
1693a40ed3dSBarry Smith   PetscFunctionReturn(0);
170b97cf49bSBarry Smith }
171b97cf49bSBarry Smith 
1724a2ae208SSatish Balay #undef __FUNCT__
1734a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
174ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
175b97cf49bSBarry Smith {
1763eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
177dfbe8321SBarry Smith   PetscErrorCode ierr;
178ace3abfcSBarry Smith   PetscBool      flag;
179b97cf49bSBarry Smith 
180433994e6SBarry Smith   PetscFunctionBegin;
181b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
182d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1830f5bd95cSBarry Smith     flag = PETSC_FALSE;
184b97cf49bSBarry Smith   }
185b97cf49bSBarry Smith 
186b97cf49bSBarry Smith   /* if the a->i are the same */
187d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
188b97cf49bSBarry Smith 
189b97cf49bSBarry Smith   /* if a->j are the same */
190b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
191b97cf49bSBarry Smith 
192ce94432eSBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1933a40ed3dSBarry Smith   PetscFunctionReturn(0);
194b97cf49bSBarry Smith }
195b97cf49bSBarry Smith 
1964a2ae208SSatish Balay #undef __FUNCT__
1974a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
1981a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
1996cd91f64SBarry Smith {
200b24ad042SBarry Smith   PetscInt   i;
2016cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2021a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
2036cd91f64SBarry Smith 
2046cd91f64SBarry Smith   PetscFunctionBegin;
205d0f46423SBarry Smith   *m    = A->rmap->n;
2066cd91f64SBarry Smith   *ia   = a->i;
2076cd91f64SBarry Smith   *ja   = a->j;
2086cd91f64SBarry Smith   *done = PETSC_TRUE;
209d42735a1SBarry Smith   if (oshift) {
210d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
211d42735a1SBarry Smith       (*ja)[i]++;
212d42735a1SBarry Smith     }
213d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
214d42735a1SBarry Smith   }
215d42735a1SBarry Smith   PetscFunctionReturn(0);
216d42735a1SBarry Smith }
217d42735a1SBarry Smith 
2184a2ae208SSatish Balay #undef __FUNCT__
2194a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
2201a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
221d42735a1SBarry Smith {
222b24ad042SBarry Smith   PetscInt   i;
223d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2241a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
225d42735a1SBarry Smith 
226d42735a1SBarry Smith   PetscFunctionBegin;
227e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
228e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
229d42735a1SBarry Smith   if (oshift) {
230d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
231d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
232d42735a1SBarry Smith       (*ja)[i]--;
233d42735a1SBarry Smith     }
234d42735a1SBarry Smith   }
2356cd91f64SBarry Smith   PetscFunctionReturn(0);
2366cd91f64SBarry Smith }
237b97cf49bSBarry Smith 
23817667f90SBarry Smith #undef __FUNCT__
23917667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
24019fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
24117667f90SBarry Smith {
24217667f90SBarry Smith   Mat               B;
24317667f90SBarry Smith   PetscErrorCode    ierr;
24417667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
24517667f90SBarry Smith   const PetscInt    *rj;
24617667f90SBarry Smith   const PetscScalar *ra;
24717667f90SBarry Smith   MPI_Comm          comm;
24817667f90SBarry Smith 
24917667f90SBarry Smith   PetscFunctionBegin;
2500298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
2510298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
2520298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
25317667f90SBarry Smith 
25417667f90SBarry Smith   /* count the number of nonzeros per row */
25517667f90SBarry Smith   for (i=0; i<m; i++) {
2560298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
2575ee9ba1cSJed Brown     for (j=0; j<len; j++) {
2585ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
2595ee9ba1cSJed Brown     }
26017667f90SBarry Smith     nzeros += len;
2610f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
26217667f90SBarry Smith   }
26317667f90SBarry Smith 
26417667f90SBarry Smith   /* malloc space for nonzeros */
265854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr);
266854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr);
267854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr);
26817667f90SBarry Smith 
26917667f90SBarry Smith   nzeros = 0;
27017667f90SBarry Smith   ia[0]  = 0;
27117667f90SBarry Smith   for (i=0; i<m; i++) {
27217667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27317667f90SBarry Smith     cnt  = 0;
27417667f90SBarry Smith     for (j=0; j<len; j++) {
2755ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
27617667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27717667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
27817667f90SBarry Smith       }
2795ee9ba1cSJed Brown     }
28017667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
28117667f90SBarry Smith     nzeros += cnt;
28217667f90SBarry Smith     ia[i+1] = nzeros;
28317667f90SBarry Smith   }
28417667f90SBarry Smith 
28517667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
28617667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
28758ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28817667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28917667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
29017667f90SBarry Smith 
29117667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
29217667f90SBarry Smith     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
29317667f90SBarry Smith   } else {
29417667f90SBarry Smith     *newmat = B;
29517667f90SBarry Smith   }
29617667f90SBarry Smith   PetscFunctionReturn(0);
29717667f90SBarry Smith }
29817667f90SBarry Smith 
299b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
3000b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
3013eda8832SBarry Smith                                        MatGetRow_MPIAdj,
3023eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
303b97cf49bSBarry Smith                                        0,
30497304618SKris Buschelman                                 /* 4*/ 0,
305b97cf49bSBarry Smith                                        0,
306b97cf49bSBarry Smith                                        0,
307b97cf49bSBarry Smith                                        0,
3080b3b1487SBarry Smith                                        0,
3090b3b1487SBarry Smith                                        0,
31097304618SKris Buschelman                                 /*10*/ 0,
3110b3b1487SBarry Smith                                        0,
3120b3b1487SBarry Smith                                        0,
3130b3b1487SBarry Smith                                        0,
3140b3b1487SBarry Smith                                        0,
31597304618SKris Buschelman                                 /*15*/ 0,
3163eda8832SBarry Smith                                        MatEqual_MPIAdj,
3170b3b1487SBarry Smith                                        0,
3180b3b1487SBarry Smith                                        0,
3190b3b1487SBarry Smith                                        0,
32097304618SKris Buschelman                                 /*20*/ 0,
3210b3b1487SBarry Smith                                        0,
3223eda8832SBarry Smith                                        MatSetOption_MPIAdj,
3230b3b1487SBarry Smith                                        0,
324d519adbfSMatthew Knepley                                 /*24*/ 0,
3250b3b1487SBarry Smith                                        0,
3260b3b1487SBarry Smith                                        0,
3270b3b1487SBarry Smith                                        0,
3280b3b1487SBarry Smith                                        0,
329d519adbfSMatthew Knepley                                 /*29*/ 0,
3300b3b1487SBarry Smith                                        0,
331273d9f13SBarry Smith                                        0,
332273d9f13SBarry Smith                                        0,
3330b3b1487SBarry Smith                                        0,
334d519adbfSMatthew Knepley                                 /*34*/ 0,
3350b3b1487SBarry Smith                                        0,
3360b3b1487SBarry Smith                                        0,
3370b3b1487SBarry Smith                                        0,
3380b3b1487SBarry Smith                                        0,
339d519adbfSMatthew Knepley                                 /*39*/ 0,
3405495542aSFande Kong                                        MatGetSubMatrices_MPIAdj,
3410b3b1487SBarry Smith                                        0,
3420b3b1487SBarry Smith                                        0,
3430b3b1487SBarry Smith                                        0,
344d519adbfSMatthew Knepley                                 /*44*/ 0,
3450b3b1487SBarry Smith                                        0,
3467d68702bSBarry Smith                                        MatShift_Basic,
3470b3b1487SBarry Smith                                        0,
3480b3b1487SBarry Smith                                        0,
349d519adbfSMatthew Knepley                                 /*49*/ 0,
3506cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
351d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
352b97cf49bSBarry Smith                                        0,
353b97cf49bSBarry Smith                                        0,
354d519adbfSMatthew Knepley                                 /*54*/ 0,
355b97cf49bSBarry Smith                                        0,
3560752156aSBarry Smith                                        0,
3570752156aSBarry Smith                                        0,
3580b3b1487SBarry Smith                                        0,
359d519adbfSMatthew Knepley                                 /*59*/ 0,
360b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
361b9b97703SBarry Smith                                        MatView_MPIAdj,
36217667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
36397304618SKris Buschelman                                        0,
364d519adbfSMatthew Knepley                                 /*64*/ 0,
36597304618SKris Buschelman                                        0,
36697304618SKris Buschelman                                        0,
36797304618SKris Buschelman                                        0,
36897304618SKris Buschelman                                        0,
369d519adbfSMatthew Knepley                                 /*69*/ 0,
37097304618SKris Buschelman                                        0,
37197304618SKris Buschelman                                        0,
37297304618SKris Buschelman                                        0,
37397304618SKris Buschelman                                        0,
374d519adbfSMatthew Knepley                                 /*74*/ 0,
37597304618SKris Buschelman                                        0,
37697304618SKris Buschelman                                        0,
37797304618SKris Buschelman                                        0,
37897304618SKris Buschelman                                        0,
379d519adbfSMatthew Knepley                                 /*79*/ 0,
38097304618SKris Buschelman                                        0,
38197304618SKris Buschelman                                        0,
38297304618SKris Buschelman                                        0,
383865e5f61SKris Buschelman                                        0,
384d519adbfSMatthew Knepley                                 /*84*/ 0,
385865e5f61SKris Buschelman                                        0,
386865e5f61SKris Buschelman                                        0,
387865e5f61SKris Buschelman                                        0,
388865e5f61SKris Buschelman                                        0,
389d519adbfSMatthew Knepley                                 /*89*/ 0,
390865e5f61SKris Buschelman                                        0,
391865e5f61SKris Buschelman                                        0,
392865e5f61SKris Buschelman                                        0,
393865e5f61SKris Buschelman                                        0,
394d519adbfSMatthew Knepley                                 /*94*/ 0,
395865e5f61SKris Buschelman                                        0,
396865e5f61SKris Buschelman                                        0,
3973964eb88SJed Brown                                        0,
3983964eb88SJed Brown                                        0,
3993964eb88SJed Brown                                 /*99*/ 0,
4003964eb88SJed Brown                                        0,
4013964eb88SJed Brown                                        0,
4023964eb88SJed Brown                                        0,
4033964eb88SJed Brown                                        0,
4043964eb88SJed Brown                                /*104*/ 0,
4053964eb88SJed Brown                                        0,
4063964eb88SJed Brown                                        0,
4073964eb88SJed Brown                                        0,
4083964eb88SJed Brown                                        0,
4093964eb88SJed Brown                                /*109*/ 0,
4103964eb88SJed Brown                                        0,
4113964eb88SJed Brown                                        0,
4123964eb88SJed Brown                                        0,
4133964eb88SJed Brown                                        0,
4143964eb88SJed Brown                                /*114*/ 0,
4153964eb88SJed Brown                                        0,
4163964eb88SJed Brown                                        0,
4173964eb88SJed Brown                                        0,
4183964eb88SJed Brown                                        0,
4193964eb88SJed Brown                                /*119*/ 0,
4203964eb88SJed Brown                                        0,
4213964eb88SJed Brown                                        0,
4223964eb88SJed Brown                                        0,
4233964eb88SJed Brown                                        0,
4243964eb88SJed Brown                                /*124*/ 0,
4253964eb88SJed Brown                                        0,
4263964eb88SJed Brown                                        0,
4273964eb88SJed Brown                                        0,
4283964eb88SJed Brown                                        0,
4293964eb88SJed Brown                                /*129*/ 0,
4303964eb88SJed Brown                                        0,
4313964eb88SJed Brown                                        0,
4323964eb88SJed Brown                                        0,
4333964eb88SJed Brown                                        0,
4343964eb88SJed Brown                                /*134*/ 0,
4353964eb88SJed Brown                                        0,
4363964eb88SJed Brown                                        0,
4373964eb88SJed Brown                                        0,
4383964eb88SJed Brown                                        0,
4393964eb88SJed Brown                                /*139*/ 0,
440f9426fe0SMark Adams                                        0,
4413964eb88SJed Brown                                        0
4423964eb88SJed Brown };
443b97cf49bSBarry Smith 
444a23d5eceSKris Buschelman #undef __FUNCT__
445a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
446f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
447a23d5eceSKris Buschelman {
448a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
449dfbe8321SBarry Smith   PetscErrorCode ierr;
4502515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
451b24ad042SBarry Smith   PetscInt ii;
452a23d5eceSKris Buschelman #endif
453a23d5eceSKris Buschelman 
454a23d5eceSKris Buschelman   PetscFunctionBegin;
45526283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
45626283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
45758ec18a5SLisandro Dalcin 
4582515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
459e32f2f54SBarry Smith   if (i[0] != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]);
460d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
461e02043d6SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]);
462a23d5eceSKris Buschelman   }
463d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
464e02043d6SBarry Smith     if (j[ii] < 0 || j[ii] >= B->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]);
465a23d5eceSKris Buschelman   }
466a23d5eceSKris Buschelman #endif
46758ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
468a23d5eceSKris Buschelman 
469a23d5eceSKris Buschelman   b->j      = j;
470a23d5eceSKris Buschelman   b->i      = i;
471a23d5eceSKris Buschelman   b->values = values;
472a23d5eceSKris Buschelman 
473d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
474a23d5eceSKris Buschelman   b->diag      = 0;
475a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
476a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
477a23d5eceSKris Buschelman 
478a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
479a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
480a23d5eceSKris Buschelman   PetscFunctionReturn(0);
481a23d5eceSKris Buschelman }
482b97cf49bSBarry Smith 
4835495542aSFande Kong static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS irows, IS icols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values);
48485374dbcSFande Kong 
4859a3665c6SJed Brown #undef __FUNCT__
48685374dbcSFande Kong #define __FUNCT__ "MatGetSubMatrices_MPIAdj"
48785374dbcSFande Kong PetscErrorCode MatGetSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
48885374dbcSFande Kong {
4895495542aSFande Kong   PetscInt           i,irow_n,icol_n,*sxadj,*sadjncy,*svalues;
4905495542aSFande Kong   PetscInt          *indices,nindx,j,k,loc;
49185374dbcSFande Kong   const PetscInt    *irow_indices,*icol_indices;
49285374dbcSFande Kong   PetscErrorCode     ierr;
49385374dbcSFande Kong 
49485374dbcSFande Kong   PetscFunctionBegin;
4955495542aSFande Kong   nindx = 0;
49685374dbcSFande Kong   for(i=0; i<n; i++){
4975495542aSFande Kong     ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr);
4985495542aSFande Kong     ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr);
4995495542aSFande Kong     nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n);
50085374dbcSFande Kong   }
5015495542aSFande Kong   ierr = PetscCalloc1(nindx,&indices);CHKERRQ(ierr);
50285374dbcSFande Kong   for(i=0; i<n; i++){
503*1b536640SFande Kong 	sxadj=0; sadjncy=0; svalues=0;
50485374dbcSFande Kong     ierr = MatGetSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr);
5055495542aSFande Kong     ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr);
5065495542aSFande Kong     ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr);
50785374dbcSFande Kong     ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr);
5085495542aSFande Kong     ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_n);CHKERRQ(ierr);
50985374dbcSFande Kong     ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr);
51085374dbcSFande Kong     ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr);
5115495542aSFande Kong     ierr = PetscMemcpy(indices+irow_n,icol_indices,sizeof(PetscInt)*icol_n);CHKERRQ(ierr);
51285374dbcSFande Kong     ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr);
5135495542aSFande Kong     nindx = irow_n+icol_n;
5145495542aSFande Kong     ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr);
515a1f6d831SFande Kong     /* renumber columns */
5165495542aSFande Kong     for(j=0; j<irow_n; j++){
517a1f6d831SFande Kong       for(k=sxadj[j]; k<sxadj[j+1]; k++){
5185495542aSFande Kong     	ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr);
51985374dbcSFande Kong #if PETSC_USE_DEBUG
52085374dbcSFande Kong     	if(loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %d \n",sadjncy[k]);
52185374dbcSFande Kong #endif
52285374dbcSFande Kong         sadjncy[k] = loc;
52385374dbcSFande Kong       }
52485374dbcSFande Kong     }
52585374dbcSFande Kong     if(scall==MAT_INITIAL_MATRIX){
5265495542aSFande Kong       ierr = MatCreateMPIAdj(PETSC_COMM_SELF,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr);
52785374dbcSFande Kong     }else{
52885374dbcSFande Kong        Mat                sadj = *(submat[i]);
52985374dbcSFande Kong        Mat_MPIAdj        *sa = (Mat_MPIAdj*)((sadj)->data);
5305495542aSFande Kong        ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_n+1));CHKERRQ(ierr);
5315495542aSFande Kong        ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);
5325495542aSFande Kong        if(svalues){ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);}
53385374dbcSFande Kong        ierr = PetscFree(sxadj);CHKERRQ(ierr);
53485374dbcSFande Kong        ierr = PetscFree(sadjncy);CHKERRQ(ierr);
5355495542aSFande Kong        if(svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);}
53685374dbcSFande Kong     }
53785374dbcSFande Kong   }
5385495542aSFande Kong   ierr = PetscFree(indices);CHKERRQ(ierr);
53985374dbcSFande Kong   PetscFunctionReturn(0);
54085374dbcSFande Kong }
54185374dbcSFande Kong 
54285374dbcSFande Kong 
54385374dbcSFande Kong #undef __FUNCT__
54485374dbcSFande Kong #define __FUNCT__ "MatGetSubMatrix_MPIAdj_data"
5455495542aSFande Kong /*
5465495542aSFande Kong  * The interface should be easy to use for both MatGetSubMatrix (parallel sub-matrix) and MatGetSubMatrices (sequential sub-matrices)
5475495542aSFande Kong  * */
5485495542aSFande Kong static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS irows, IS icols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values)
549841d17a1SFande Kong {
5505495542aSFande Kong   PetscInt        	 nlrows_is,icols_n,i,j,nroots,nleaves,owner,rlocalindex,*ncols_send,*ncols_recv;
5515495542aSFande Kong   PetscInt           nlrows_mat,*adjncy_recv,Ncols_recv,Ncols_send,*xadj_recv,*values_recv;
5525495542aSFande Kong   PetscInt          *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues,isvalue;
5535495542aSFande Kong   const PetscInt    *irows_indices,*icols_indices,*xadj, *adjncy;
554841d17a1SFande Kong   Mat_MPIAdj        *a = (Mat_MPIAdj*)adj->data;
555841d17a1SFande Kong   PetscLayout        rmap;
556841d17a1SFande Kong   MPI_Comm           comm;
557841d17a1SFande Kong   PetscSF            sf;
558841d17a1SFande Kong   PetscSFNode       *iremote;
559841d17a1SFande Kong   PetscBool          done;
560841d17a1SFande Kong   PetscErrorCode     ierr;
561841d17a1SFande Kong 
562841d17a1SFande Kong   PetscFunctionBegin;
563841d17a1SFande Kong   /* communicator */
564841d17a1SFande Kong   ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr);
565841d17a1SFande Kong   /* Layouts */
566841d17a1SFande Kong   ierr = MatGetLayouts(adj,&rmap,PETSC_NULL);CHKERRQ(ierr);
567841d17a1SFande Kong   /* get rows information */
5685495542aSFande Kong   ierr = ISGetLocalSize(irows,&nlrows_is);CHKERRQ(ierr);
5695495542aSFande Kong   ierr = ISGetIndices(irows,&irows_indices);CHKERRQ(ierr);
5705495542aSFande Kong   ierr = PetscCalloc1(nlrows_is,&iremote);CHKERRQ(ierr);
5715495542aSFande Kong   /* construct sf graph*/
5725495542aSFande Kong   nleaves = nlrows_is;
5735495542aSFande Kong   for(i=0; i<nlrows_is; i++){
5742253a2e4SFande Kong 	owner = -1;
5752253a2e4SFande Kong 	rlocalindex = -1;
5765495542aSFande Kong     ierr = PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex);CHKERRQ(ierr);
577841d17a1SFande Kong     iremote[i].rank  = owner;
5785495542aSFande Kong     iremote[i].index = rlocalindex;
579841d17a1SFande Kong   }
5805495542aSFande Kong   ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr);
5815495542aSFande Kong   ierr = PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv);CHKERRQ(ierr);
5825495542aSFande Kong   nroots = nlrows_mat;
5835495542aSFande Kong   for(i=0; i<nlrows_mat; i++){
584841d17a1SFande Kong 	ncols_send[i] = xadj[i+1]-xadj[i];
585841d17a1SFande Kong   }
586841d17a1SFande Kong   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
58785374dbcSFande Kong   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
588841d17a1SFande Kong   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
589841d17a1SFande Kong   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
590841d17a1SFande Kong   ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr);
591841d17a1SFande Kong   ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr);
592841d17a1SFande Kong   ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr);
593841d17a1SFande Kong   ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr);
594841d17a1SFande Kong   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
595841d17a1SFande Kong   Ncols_recv =0;
5965495542aSFande Kong   for(i=0; i<nlrows_is; i++){
597841d17a1SFande Kong 	 Ncols_recv             += ncols_recv[i];
59885374dbcSFande Kong 	 ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i];
599841d17a1SFande Kong   }
600841d17a1SFande Kong   Ncols_send = 0;
6015495542aSFande Kong   for(i=0; i<nlrows_mat; i++){
602841d17a1SFande Kong 	Ncols_send += ncols_send[i];
603841d17a1SFande Kong   }
604841d17a1SFande Kong   ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr);
6055495542aSFande Kong   ierr = PetscCalloc1(Ncols_recv,&adjncy_recv);CHKERRQ(ierr);
606841d17a1SFande Kong   nleaves = Ncols_recv;
6075495542aSFande Kong   Ncols_recv = 0;
6085495542aSFande Kong   for(i=0; i<nlrows_is; i++){
6095495542aSFande Kong     ierr = PetscLayoutFindOwner(rmap,irows_indices[i],&owner);CHKERRQ(ierr);
610841d17a1SFande Kong     for(j=0; j<ncols_recv[i]; j++){
6115495542aSFande Kong       iremote[Ncols_recv].rank    = owner;
6125495542aSFande Kong       iremote[Ncols_recv++].index = xadj_recv[i]+j;
613841d17a1SFande Kong     }
614841d17a1SFande Kong   }
6155495542aSFande Kong   ierr = ISRestoreIndices(irows,&irows_indices);CHKERRQ(ierr);
616a1f6d831SFande Kong   /*if we need to deal with edge weights ???*/
6175495542aSFande Kong   if(a->values){isvalue=1;}else{isvalue=0;}
6185495542aSFande Kong   /*involve a global communication */
619a1f6d831SFande Kong   /*ierr = MPI_Allreduce(&isvalue,&isvalue,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);*/
6205495542aSFande Kong   if(isvalue){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);}
621841d17a1SFande Kong   nroots = Ncols_send;
622841d17a1SFande Kong   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
62385374dbcSFande Kong   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
624841d17a1SFande Kong   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
625841d17a1SFande Kong   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
626841d17a1SFande Kong   ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr);
627841d17a1SFande Kong   ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr);
6285495542aSFande Kong   if(isvalue){
629841d17a1SFande Kong 	ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr);
630841d17a1SFande Kong 	ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr);
6315495542aSFande Kong   }
632841d17a1SFande Kong   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
6335495542aSFande Kong   ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr);
6345495542aSFande Kong   ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr);
6355495542aSFande Kong   ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr);
63685374dbcSFande Kong   rnclos = 0;
6375495542aSFande Kong   for(i=0; i<nlrows_is; i++){
63885374dbcSFande Kong     for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){
6395495542aSFande Kong       ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr);
64085374dbcSFande Kong       if(loc<0){
64185374dbcSFande Kong         adjncy_recv[j] = -1;
6425495542aSFande Kong         if(isvalue) values_recv[j] = -1;
64385374dbcSFande Kong         ncols_recv[i]--;
64485374dbcSFande Kong       }else{
64585374dbcSFande Kong     	rnclos++;
64685374dbcSFande Kong       }
64785374dbcSFande Kong     }
64885374dbcSFande Kong   }
6495495542aSFande Kong   ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr);
65085374dbcSFande Kong   ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr);
6515495542aSFande Kong   if(isvalue) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);}
6525495542aSFande Kong   ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr);
65385374dbcSFande Kong   rnclos = 0;
6545495542aSFande Kong   for(i=0; i<nlrows_is; i++){
65585374dbcSFande Kong 	for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){
65685374dbcSFande Kong 	  if(adjncy_recv[j]<0) continue;
65785374dbcSFande Kong 	  sadjncy[rnclos] = adjncy_recv[j];
6585495542aSFande Kong 	  if(isvalue) svalues[rnclos] = values_recv[j];
65985374dbcSFande Kong 	  rnclos++;
66085374dbcSFande Kong 	}
66185374dbcSFande Kong   }
6625495542aSFande Kong   for(i=0; i<nlrows_is; i++){
66385374dbcSFande Kong 	sxadj[i+1] = sxadj[i]+ncols_recv[i];
66485374dbcSFande Kong   }
66585374dbcSFande Kong   if(sadj_xadj)  { *sadj_xadj = sxadj;}else    { ierr = PetscFree(sxadj);CHKERRQ(ierr);}
66685374dbcSFande Kong   if(sadj_adjncy){ *sadj_adjncy = sadjncy;}else{ ierr = PetscFree(sadjncy);CHKERRQ(ierr);}
6675495542aSFande Kong   if(sadj_values){
6685495542aSFande Kong 	if(isvalue) *sadj_values = svalues; else *sadj_values=0;
6695495542aSFande Kong   }else{
6705495542aSFande Kong 	if(isvalue) {ierr = PetscFree(svalues);CHKERRQ(ierr);}
6715495542aSFande Kong   }
67285374dbcSFande Kong   ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr);
6735495542aSFande Kong   ierr = PetscFree(adjncy_recv);CHKERRQ(ierr);
6745495542aSFande Kong   if(isvalue) {ierr = PetscFree(values_recv);CHKERRQ(ierr);}
675841d17a1SFande Kong   PetscFunctionReturn(0);
676841d17a1SFande Kong }
677841d17a1SFande Kong 
678841d17a1SFande Kong 
679841d17a1SFande Kong #undef __FUNCT__
6809a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj"
681f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
6829a3665c6SJed Brown {
6839a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
6849a3665c6SJed Brown   PetscErrorCode ierr;
6859a3665c6SJed Brown   const PetscInt *ranges;
6869a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
6879a3665c6SJed Brown   MPI_Group      agroup,bgroup;
6889a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
6899a3665c6SJed Brown 
6909a3665c6SJed Brown   PetscFunctionBegin;
6910298fd71SBarry Smith   *B    = NULL;
692ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
6939a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
6949a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
6959a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
6969a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
6979a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
6989a3665c6SJed Brown   }
6999a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
7009a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
7019a3665c6SJed Brown     *B   = A;
7029a3665c6SJed Brown     PetscFunctionReturn(0);
7039a3665c6SJed Brown   }
7049a3665c6SJed Brown 
705785e854fSJed Brown   ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr);
7069a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
7079a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
7089a3665c6SJed Brown   }
7099a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
7109a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
7119a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
7129a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
7139a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
7149a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
7159a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
7169a3665c6SJed Brown     PetscInt   m,N;
7179a3665c6SJed Brown     Mat_MPIAdj *b;
7180298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
7190298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
7209a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
7219a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
7229a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
7239a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
7249a3665c6SJed Brown   }
7259a3665c6SJed Brown   PetscFunctionReturn(0);
7269a3665c6SJed Brown }
7279a3665c6SJed Brown 
7289a3665c6SJed Brown #undef __FUNCT__
7299a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat"
7309a3665c6SJed Brown /*@
7319a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
7329a3665c6SJed Brown 
7339a3665c6SJed Brown    Collective
7349a3665c6SJed Brown 
7359a3665c6SJed Brown    Input Arguments:
7369a3665c6SJed Brown .  A - original MPIAdj matrix
7379a3665c6SJed Brown 
7389a3665c6SJed Brown    Output Arguments:
7390298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
7409a3665c6SJed Brown 
7419a3665c6SJed Brown    Level: developer
7429a3665c6SJed Brown 
7439a3665c6SJed Brown    Note:
7449a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
7459a3665c6SJed Brown 
7469a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
7479a3665c6SJed Brown 
7489a3665c6SJed Brown .seealso: MatCreateMPIAdj()
7499a3665c6SJed Brown @*/
7509a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
7519a3665c6SJed Brown {
7529a3665c6SJed Brown   PetscErrorCode ierr;
7539a3665c6SJed Brown 
7549a3665c6SJed Brown   PetscFunctionBegin;
7559a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
7569a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
7579a3665c6SJed Brown   PetscFunctionReturn(0);
7589a3665c6SJed Brown }
7599a3665c6SJed Brown 
7600bad9183SKris Buschelman /*MC
761fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
7620bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
7630bad9183SKris Buschelman 
7640bad9183SKris Buschelman   Level: beginner
7650bad9183SKris Buschelman 
7660bad9183SKris Buschelman .seealso: MatCreateMPIAdj
7670bad9183SKris Buschelman M*/
7680bad9183SKris Buschelman 
7694a2ae208SSatish Balay #undef __FUNCT__
7704a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
7718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
772273d9f13SBarry Smith {
773273d9f13SBarry Smith   Mat_MPIAdj     *b;
7746849ba73SBarry Smith   PetscErrorCode ierr;
775273d9f13SBarry Smith 
776273d9f13SBarry Smith   PetscFunctionBegin;
777b00a9115SJed Brown   ierr         = PetscNewLog(B,&b);CHKERRQ(ierr);
778b0a32e0cSBarry Smith   B->data      = (void*)b;
779273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
780273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
781273d9f13SBarry Smith 
782bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
783bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
78417667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
785273d9f13SBarry Smith   PetscFunctionReturn(0);
786273d9f13SBarry Smith }
787273d9f13SBarry Smith 
7884a2ae208SSatish Balay #undef __FUNCT__
7894a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
790a23d5eceSKris Buschelman /*@C
791a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
792a23d5eceSKris Buschelman 
7933f9fe445SBarry Smith    Logically Collective on MPI_Comm
794a23d5eceSKris Buschelman 
795a23d5eceSKris Buschelman    Input Parameters:
796a23d5eceSKris Buschelman +  A - the matrix
797a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
798a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
799a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
800a23d5eceSKris Buschelman -  values - [optional] edge weights
801a23d5eceSKris Buschelman 
802a23d5eceSKris Buschelman    Level: intermediate
803a23d5eceSKris Buschelman 
804a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
805a23d5eceSKris Buschelman @*/
8067087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
807273d9f13SBarry Smith {
8084ac538c5SBarry Smith   PetscErrorCode ierr;
809273d9f13SBarry Smith 
810273d9f13SBarry Smith   PetscFunctionBegin;
8114ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
812273d9f13SBarry Smith   PetscFunctionReturn(0);
813273d9f13SBarry Smith }
814273d9f13SBarry Smith 
8154a2ae208SSatish Balay #undef __FUNCT__
8164a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
817b97cf49bSBarry Smith /*@C
8183eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
819b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
820b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
821b97cf49bSBarry Smith 
822ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
823ef5ee4d1SLois Curfman McInnes 
824b97cf49bSBarry Smith    Input Parameters:
825c2e958feSBarry Smith +  comm - MPI communicator
8260752156aSBarry Smith .  m - number of local rows
82758ec18a5SLisandro Dalcin .  N - number of global columns
828b97cf49bSBarry Smith .  i - the indices into j for the start of each row
829329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
830ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
831329f5518SBarry Smith -  values -[optional] edge weights
832b97cf49bSBarry Smith 
833b97cf49bSBarry Smith    Output Parameter:
834b97cf49bSBarry Smith .  A - the matrix
835b97cf49bSBarry Smith 
83615091d37SBarry Smith    Level: intermediate
83715091d37SBarry Smith 
8384bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
8394bc6d8c0SBarry Smith    MatSetValues().
840c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
841ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
8421198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
843327e1a73SBarry Smith    Should not include the matrix diagonals.
844b97cf49bSBarry Smith 
845e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
846ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
847ededeb1aSvictorle 
848ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
849b97cf49bSBarry Smith 
850e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
851b97cf49bSBarry Smith @*/
8527087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
853b97cf49bSBarry Smith {
854dfbe8321SBarry Smith   PetscErrorCode ierr;
855b97cf49bSBarry Smith 
856433994e6SBarry Smith   PetscFunctionBegin;
857f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
85858ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
859273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
860273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
8613a40ed3dSBarry Smith   PetscFunctionReturn(0);
862b97cf49bSBarry Smith }
863b97cf49bSBarry Smith 
864c2e958feSBarry Smith 
865433994e6SBarry Smith 
866433994e6SBarry Smith 
867433994e6SBarry Smith 
868433994e6SBarry Smith 
869433994e6SBarry Smith 
870