xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision b2566f29c2b6470df769aa9f7deb9e2726b0959e)
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);
236c4ed002SBarry Smith   } else if (format == PETSC_VIEWER_ASCII_MATLAB) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported");
246c4ed002SBarry Smith   else {
25d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
26c5e4d11fSDmitry Karpeev     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
27b97cf49bSBarry Smith     for (i=0; i<m; i++) {
28d0f46423SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr);
29b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
3077431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
310752156aSBarry Smith       }
32b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
33b97cf49bSBarry Smith     }
34d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
35b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
36c5e4d11fSDmitry Karpeev     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
377b23a99aSBarry Smith   }
383a40ed3dSBarry Smith   PetscFunctionReturn(0);
39b97cf49bSBarry Smith }
40b97cf49bSBarry Smith 
414a2ae208SSatish Balay #undef __FUNCT__
424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
43dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
44b97cf49bSBarry Smith {
45dfbe8321SBarry Smith   PetscErrorCode ierr;
46ace3abfcSBarry Smith   PetscBool      iascii;
47b97cf49bSBarry Smith 
48433994e6SBarry Smith   PetscFunctionBegin;
49251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5032077d6dSBarry Smith   if (iascii) {
513eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
52b97cf49bSBarry Smith   }
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54b97cf49bSBarry Smith }
55b97cf49bSBarry Smith 
564a2ae208SSatish Balay #undef __FUNCT__
574a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
58dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat)
59b97cf49bSBarry Smith {
603eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
61dfbe8321SBarry Smith   PetscErrorCode ierr;
6294d884c6SBarry Smith 
6394d884c6SBarry Smith   PetscFunctionBegin;
64aa482453SBarry Smith #if defined(PETSC_USE_LOG)
65d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
66b97cf49bSBarry Smith #endif
6705b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
680400557aSBarry Smith   if (a->freeaij) {
6914196391SBarry Smith     if (a->freeaijwithfree) {
7014196391SBarry Smith       if (a->i) free(a->i);
7114196391SBarry Smith       if (a->j) free(a->j);
7214196391SBarry Smith     } else {
73606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
74606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
7505b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
7614196391SBarry Smith     }
770400557aSBarry Smith   }
782b1d8763SJed Brown   ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
79bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
80dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
81bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
82bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
833a40ed3dSBarry Smith   PetscFunctionReturn(0);
84b97cf49bSBarry Smith }
85b97cf49bSBarry Smith 
864a2ae208SSatish Balay #undef __FUNCT__
874a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
88ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
89b97cf49bSBarry Smith {
903eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
9163ba0a88SBarry Smith   PetscErrorCode ierr;
92b97cf49bSBarry Smith 
93433994e6SBarry Smith   PetscFunctionBegin;
9412c028f9SKris Buschelman   switch (op) {
9577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9612c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
979a4540c5SBarry Smith   case MAT_HERMITIAN:
984e0d8c25SBarry Smith     a->symmetric = flg;
999a4540c5SBarry Smith     break;
1009a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1019a4540c5SBarry Smith     break;
10212c028f9SKris Buschelman   default:
103290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
10412c028f9SKris Buschelman     break;
105b97cf49bSBarry Smith   }
1063a40ed3dSBarry Smith   PetscFunctionReturn(0);
107b97cf49bSBarry Smith }
108b97cf49bSBarry Smith 
109b97cf49bSBarry Smith 
110b97cf49bSBarry Smith /*
111b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
112b97cf49bSBarry Smith */
113b97cf49bSBarry Smith 
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
116dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
117b97cf49bSBarry Smith {
1183eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1196849ba73SBarry Smith   PetscErrorCode ierr;
120d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n;
121b97cf49bSBarry Smith 
122433994e6SBarry Smith   PetscFunctionBegin;
123785e854fSJed Brown   ierr = PetscMalloc1(m,&a->diag);CHKERRQ(ierr);
1243bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,m*sizeof(PetscInt));CHKERRQ(ierr);
125d0f46423SBarry Smith   for (i=0; i<A->rmap->n; i++) {
126b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
127b97cf49bSBarry Smith       if (a->j[j] == i) {
12809f38230SBarry Smith         a->diag[i] = j;
129b97cf49bSBarry Smith         break;
130b97cf49bSBarry Smith       }
131b97cf49bSBarry Smith     }
132b97cf49bSBarry Smith   }
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
134b97cf49bSBarry Smith }
135b97cf49bSBarry Smith 
1364a2ae208SSatish Balay #undef __FUNCT__
1374a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
138b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
139b97cf49bSBarry Smith {
1403eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
1412b1d8763SJed Brown   PetscErrorCode ierr;
142b97cf49bSBarry Smith 
143433994e6SBarry Smith   PetscFunctionBegin;
144d0f46423SBarry Smith   row -= A->rmap->rstart;
145e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
146b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
1472b1d8763SJed Brown   if (v) {
1482b1d8763SJed Brown     PetscInt j;
1492b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
1502b1d8763SJed Brown       ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
1512b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz);
152785e854fSJed Brown       ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr);
153b97cf49bSBarry Smith     }
15491f8cafeSFande Kong     for (j=0; j<*nz; j++){
15591f8cafeSFande Kong       a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0;
15691f8cafeSFande Kong     }
1572b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
158b97cf49bSBarry Smith   }
15992b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161b97cf49bSBarry Smith }
162b97cf49bSBarry Smith 
1634a2ae208SSatish Balay #undef __FUNCT__
1644a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
165b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
166b97cf49bSBarry Smith {
167433994e6SBarry Smith   PetscFunctionBegin;
1683a40ed3dSBarry Smith   PetscFunctionReturn(0);
169b97cf49bSBarry Smith }
170b97cf49bSBarry Smith 
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
173ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
174b97cf49bSBarry Smith {
1753eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
176dfbe8321SBarry Smith   PetscErrorCode ierr;
177ace3abfcSBarry Smith   PetscBool      flag;
178b97cf49bSBarry Smith 
179433994e6SBarry Smith   PetscFunctionBegin;
180b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
181d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1820f5bd95cSBarry Smith     flag = PETSC_FALSE;
183b97cf49bSBarry Smith   }
184b97cf49bSBarry Smith 
185b97cf49bSBarry Smith   /* if the a->i are the same */
186d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
187b97cf49bSBarry Smith 
188b97cf49bSBarry Smith   /* if a->j are the same */
189b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
190b97cf49bSBarry Smith 
191*b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
193b97cf49bSBarry Smith }
194b97cf49bSBarry Smith 
1954a2ae208SSatish Balay #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
1971a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
1986cd91f64SBarry Smith {
199b24ad042SBarry Smith   PetscInt   i;
2006cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2011a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
2026cd91f64SBarry Smith 
2036cd91f64SBarry Smith   PetscFunctionBegin;
204d0f46423SBarry Smith   *m    = A->rmap->n;
2056cd91f64SBarry Smith   *ia   = a->i;
2066cd91f64SBarry Smith   *ja   = a->j;
2076cd91f64SBarry Smith   *done = PETSC_TRUE;
208d42735a1SBarry Smith   if (oshift) {
209d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
210d42735a1SBarry Smith       (*ja)[i]++;
211d42735a1SBarry Smith     }
212d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
213d42735a1SBarry Smith   }
214d42735a1SBarry Smith   PetscFunctionReturn(0);
215d42735a1SBarry Smith }
216d42735a1SBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
2191a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
220d42735a1SBarry Smith {
221b24ad042SBarry Smith   PetscInt   i;
222d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2231a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
224d42735a1SBarry Smith 
225d42735a1SBarry Smith   PetscFunctionBegin;
226e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
227e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
228d42735a1SBarry Smith   if (oshift) {
229d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
230d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
231d42735a1SBarry Smith       (*ja)[i]--;
232d42735a1SBarry Smith     }
233d42735a1SBarry Smith   }
2346cd91f64SBarry Smith   PetscFunctionReturn(0);
2356cd91f64SBarry Smith }
236b97cf49bSBarry Smith 
23717667f90SBarry Smith #undef __FUNCT__
23817667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
23919fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
24017667f90SBarry Smith {
24117667f90SBarry Smith   Mat               B;
24217667f90SBarry Smith   PetscErrorCode    ierr;
24317667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
24417667f90SBarry Smith   const PetscInt    *rj;
24517667f90SBarry Smith   const PetscScalar *ra;
24617667f90SBarry Smith   MPI_Comm          comm;
24717667f90SBarry Smith 
24817667f90SBarry Smith   PetscFunctionBegin;
2490298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
2500298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
2510298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
25217667f90SBarry Smith 
25317667f90SBarry Smith   /* count the number of nonzeros per row */
25417667f90SBarry Smith   for (i=0; i<m; i++) {
2550298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
2565ee9ba1cSJed Brown     for (j=0; j<len; j++) {
2575ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
2585ee9ba1cSJed Brown     }
25917667f90SBarry Smith     nzeros += len;
2600f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
26117667f90SBarry Smith   }
26217667f90SBarry Smith 
26317667f90SBarry Smith   /* malloc space for nonzeros */
264854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr);
265854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr);
266854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr);
26717667f90SBarry Smith 
26817667f90SBarry Smith   nzeros = 0;
26917667f90SBarry Smith   ia[0]  = 0;
27017667f90SBarry Smith   for (i=0; i<m; i++) {
27117667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27217667f90SBarry Smith     cnt  = 0;
27317667f90SBarry Smith     for (j=0; j<len; j++) {
2745ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
27517667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27617667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
27717667f90SBarry Smith       }
2785ee9ba1cSJed Brown     }
27917667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
28017667f90SBarry Smith     nzeros += cnt;
28117667f90SBarry Smith     ia[i+1] = nzeros;
28217667f90SBarry Smith   }
28317667f90SBarry Smith 
28417667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
28517667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
28658ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28717667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28817667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
28917667f90SBarry Smith 
29017667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
29128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
29217667f90SBarry Smith   } else {
29317667f90SBarry Smith     *newmat = B;
29417667f90SBarry Smith   }
29517667f90SBarry Smith   PetscFunctionReturn(0);
29617667f90SBarry Smith }
29717667f90SBarry Smith 
298b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2990b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
3003eda8832SBarry Smith                                        MatGetRow_MPIAdj,
3013eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
302b97cf49bSBarry Smith                                        0,
30397304618SKris Buschelman                                 /* 4*/ 0,
304b97cf49bSBarry Smith                                        0,
305b97cf49bSBarry Smith                                        0,
306b97cf49bSBarry Smith                                        0,
3070b3b1487SBarry Smith                                        0,
3080b3b1487SBarry Smith                                        0,
30997304618SKris Buschelman                                 /*10*/ 0,
3100b3b1487SBarry Smith                                        0,
3110b3b1487SBarry Smith                                        0,
3120b3b1487SBarry Smith                                        0,
3130b3b1487SBarry Smith                                        0,
31497304618SKris Buschelman                                 /*15*/ 0,
3153eda8832SBarry Smith                                        MatEqual_MPIAdj,
3160b3b1487SBarry Smith                                        0,
3170b3b1487SBarry Smith                                        0,
3180b3b1487SBarry Smith                                        0,
31997304618SKris Buschelman                                 /*20*/ 0,
3200b3b1487SBarry Smith                                        0,
3213eda8832SBarry Smith                                        MatSetOption_MPIAdj,
3220b3b1487SBarry Smith                                        0,
323d519adbfSMatthew Knepley                                 /*24*/ 0,
3240b3b1487SBarry Smith                                        0,
3250b3b1487SBarry Smith                                        0,
3260b3b1487SBarry Smith                                        0,
3270b3b1487SBarry Smith                                        0,
328d519adbfSMatthew Knepley                                 /*29*/ 0,
3290b3b1487SBarry Smith                                        0,
330273d9f13SBarry Smith                                        0,
331273d9f13SBarry Smith                                        0,
3320b3b1487SBarry Smith                                        0,
333d519adbfSMatthew Knepley                                 /*34*/ 0,
3340b3b1487SBarry Smith                                        0,
3350b3b1487SBarry Smith                                        0,
3360b3b1487SBarry Smith                                        0,
3370b3b1487SBarry Smith                                        0,
338d519adbfSMatthew Knepley                                 /*39*/ 0,
3395495542aSFande Kong                                        MatGetSubMatrices_MPIAdj,
3400b3b1487SBarry Smith                                        0,
3410b3b1487SBarry Smith                                        0,
3420b3b1487SBarry Smith                                        0,
343d519adbfSMatthew Knepley                                 /*44*/ 0,
3440b3b1487SBarry Smith                                        0,
3457d68702bSBarry Smith                                        MatShift_Basic,
3460b3b1487SBarry Smith                                        0,
3470b3b1487SBarry Smith                                        0,
348d519adbfSMatthew Knepley                                 /*49*/ 0,
3496cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
350d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
351b97cf49bSBarry Smith                                        0,
352b97cf49bSBarry Smith                                        0,
353d519adbfSMatthew Knepley                                 /*54*/ 0,
354b97cf49bSBarry Smith                                        0,
3550752156aSBarry Smith                                        0,
3560752156aSBarry Smith                                        0,
3570b3b1487SBarry Smith                                        0,
358d519adbfSMatthew Knepley                                 /*59*/ 0,
359b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
360b9b97703SBarry Smith                                        MatView_MPIAdj,
36117667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
36297304618SKris Buschelman                                        0,
363d519adbfSMatthew Knepley                                 /*64*/ 0,
36497304618SKris Buschelman                                        0,
36597304618SKris Buschelman                                        0,
36697304618SKris Buschelman                                        0,
36797304618SKris Buschelman                                        0,
368d519adbfSMatthew Knepley                                 /*69*/ 0,
36997304618SKris Buschelman                                        0,
37097304618SKris Buschelman                                        0,
37197304618SKris Buschelman                                        0,
37297304618SKris Buschelman                                        0,
373d519adbfSMatthew Knepley                                 /*74*/ 0,
37497304618SKris Buschelman                                        0,
37597304618SKris Buschelman                                        0,
37697304618SKris Buschelman                                        0,
37797304618SKris Buschelman                                        0,
378d519adbfSMatthew Knepley                                 /*79*/ 0,
37997304618SKris Buschelman                                        0,
38097304618SKris Buschelman                                        0,
38197304618SKris Buschelman                                        0,
382865e5f61SKris Buschelman                                        0,
383d519adbfSMatthew Knepley                                 /*84*/ 0,
384865e5f61SKris Buschelman                                        0,
385865e5f61SKris Buschelman                                        0,
386865e5f61SKris Buschelman                                        0,
387865e5f61SKris Buschelman                                        0,
388d519adbfSMatthew Knepley                                 /*89*/ 0,
389865e5f61SKris Buschelman                                        0,
390865e5f61SKris Buschelman                                        0,
391865e5f61SKris Buschelman                                        0,
392865e5f61SKris Buschelman                                        0,
393d519adbfSMatthew Knepley                                 /*94*/ 0,
394865e5f61SKris Buschelman                                        0,
395865e5f61SKris Buschelman                                        0,
3963964eb88SJed Brown                                        0,
3973964eb88SJed Brown                                        0,
3983964eb88SJed Brown                                 /*99*/ 0,
3993964eb88SJed Brown                                        0,
4003964eb88SJed Brown                                        0,
4013964eb88SJed Brown                                        0,
4023964eb88SJed Brown                                        0,
4033964eb88SJed Brown                                /*104*/ 0,
4043964eb88SJed Brown                                        0,
4053964eb88SJed Brown                                        0,
4063964eb88SJed Brown                                        0,
4073964eb88SJed Brown                                        0,
4083964eb88SJed Brown                                /*109*/ 0,
4093964eb88SJed Brown                                        0,
4103964eb88SJed Brown                                        0,
4113964eb88SJed Brown                                        0,
4123964eb88SJed Brown                                        0,
4133964eb88SJed Brown                                /*114*/ 0,
4143964eb88SJed Brown                                        0,
4153964eb88SJed Brown                                        0,
4163964eb88SJed Brown                                        0,
4173964eb88SJed Brown                                        0,
4183964eb88SJed Brown                                /*119*/ 0,
4193964eb88SJed Brown                                        0,
4203964eb88SJed Brown                                        0,
4213964eb88SJed Brown                                        0,
4223964eb88SJed Brown                                        0,
4233964eb88SJed Brown                                /*124*/ 0,
4243964eb88SJed Brown                                        0,
4253964eb88SJed Brown                                        0,
4263964eb88SJed Brown                                        0,
42740475a53SFande Kong                                        MatGetSubMatricesMPI_MPIAdj,
4283964eb88SJed Brown                                /*129*/ 0,
4293964eb88SJed Brown                                        0,
4303964eb88SJed Brown                                        0,
4313964eb88SJed Brown                                        0,
4323964eb88SJed Brown                                        0,
4333964eb88SJed Brown                                /*134*/ 0,
4343964eb88SJed Brown                                        0,
4353964eb88SJed Brown                                        0,
4363964eb88SJed Brown                                        0,
4373964eb88SJed Brown                                        0,
4383964eb88SJed Brown                                /*139*/ 0,
439f9426fe0SMark Adams                                        0,
4403964eb88SJed Brown                                        0
4413964eb88SJed Brown };
442b97cf49bSBarry Smith 
443a23d5eceSKris Buschelman #undef __FUNCT__
444a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
445f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
446a23d5eceSKris Buschelman {
447a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
448dfbe8321SBarry Smith   PetscErrorCode ierr;
4492515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
450b24ad042SBarry Smith   PetscInt ii;
451a23d5eceSKris Buschelman #endif
452a23d5eceSKris Buschelman 
453a23d5eceSKris Buschelman   PetscFunctionBegin;
45426283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
45526283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
45658ec18a5SLisandro Dalcin 
4572515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
458e32f2f54SBarry 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]);
459d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
460e02043d6SBarry 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]);
461a23d5eceSKris Buschelman   }
462d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
463e02043d6SBarry 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]);
464a23d5eceSKris Buschelman   }
465a23d5eceSKris Buschelman #endif
46658ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
467a23d5eceSKris Buschelman 
468a23d5eceSKris Buschelman   b->j      = j;
469a23d5eceSKris Buschelman   b->i      = i;
470a23d5eceSKris Buschelman   b->values = values;
471a23d5eceSKris Buschelman 
472d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
473a23d5eceSKris Buschelman   b->diag      = 0;
474a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
475a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
476a23d5eceSKris Buschelman 
477a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
478a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
479a23d5eceSKris Buschelman   PetscFunctionReturn(0);
480a23d5eceSKris Buschelman }
481b97cf49bSBarry Smith 
48228b21874SFande Kong static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat,IS,IS, PetscInt **,PetscInt **,PetscInt **);
48328b21874SFande Kong static PetscErrorCode MatGetSubMatrices_MPIAdj_Private(Mat,PetscInt,const IS[],const IS[],PetscBool,MatReuse,Mat **);
48440475a53SFande Kong 
48540475a53SFande Kong #undef __FUNCT__
48640475a53SFande Kong #define __FUNCT__ "MatGetSubMatricesMPI_MPIAdj"
48740475a53SFande Kong PetscErrorCode MatGetSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
48840475a53SFande Kong {
48940475a53SFande Kong   PetscErrorCode     ierr;
49040475a53SFande Kong   /*get sub-matrices across a sub communicator */
49140475a53SFande Kong   PetscFunctionBegin;
49240475a53SFande Kong   ierr = MatGetSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr);
49340475a53SFande Kong   PetscFunctionReturn(0);
49440475a53SFande Kong }
49540475a53SFande Kong 
49685374dbcSFande Kong 
4979a3665c6SJed Brown #undef __FUNCT__
49885374dbcSFande Kong #define __FUNCT__ "MatGetSubMatrices_MPIAdj"
49985374dbcSFande Kong PetscErrorCode MatGetSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
50085374dbcSFande Kong {
50140475a53SFande Kong   PetscErrorCode     ierr;
50240475a53SFande Kong 
50340475a53SFande Kong   PetscFunctionBegin;
50440475a53SFande Kong   /*get sub-matrices based on PETSC_COMM_SELF */
50540475a53SFande Kong   ierr = MatGetSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr);
50640475a53SFande Kong   PetscFunctionReturn(0);
50740475a53SFande Kong }
50840475a53SFande Kong 
50940475a53SFande Kong #undef __FUNCT__
51040475a53SFande Kong #define __FUNCT__ "MatGetSubMatrices_MPIAdj_Private"
51140475a53SFande Kong static PetscErrorCode MatGetSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[])
51240475a53SFande Kong {
5135495542aSFande Kong   PetscInt           i,irow_n,icol_n,*sxadj,*sadjncy,*svalues;
51428b21874SFande Kong   PetscInt          *indices,nindx,j,k,loc;
51528b21874SFande Kong   PetscMPIInt        issame;
51685374dbcSFande Kong   const PetscInt    *irow_indices,*icol_indices;
51740475a53SFande Kong   MPI_Comm           scomm_row,scomm_col,scomm_mat;
51885374dbcSFande Kong   PetscErrorCode     ierr;
51985374dbcSFande Kong 
52085374dbcSFande Kong   PetscFunctionBegin;
5215495542aSFande Kong   nindx = 0;
52240475a53SFande Kong   /*
52340475a53SFande Kong    * Estimate a maximum number for allocating memory
52440475a53SFande Kong    */
52585374dbcSFande Kong   for(i=0; i<n; i++){
5265495542aSFande Kong     ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr);
5275495542aSFande Kong     ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr);
5285495542aSFande Kong     nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n);
52985374dbcSFande Kong   }
5305495542aSFande Kong   ierr = PetscCalloc1(nindx,&indices);CHKERRQ(ierr);
53140475a53SFande Kong   /* construct a submat */
53285374dbcSFande Kong   for(i=0; i<n; i++){
53340475a53SFande Kong 	/*comms */
53440475a53SFande Kong     if(subcomm){
53540475a53SFande Kong 	  ierr = PetscObjectGetComm((PetscObject)irow[i],&scomm_row);CHKERRQ(ierr);
53640475a53SFande Kong 	  ierr = PetscObjectGetComm((PetscObject)icol[i],&scomm_col);CHKERRQ(ierr);
53740475a53SFande Kong 	  ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRQ(ierr);
53840475a53SFande Kong 	  if(issame != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row index set must have the same comm as the col index set\n");
53940475a53SFande Kong 	  ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRQ(ierr);
54040475a53SFande Kong 	  if(issame == MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP," can not use PETSC_COMM_SELF as comm when extracting a parallel submatrix\n");
54140475a53SFande Kong 	}else{
54240475a53SFande Kong 	  scomm_row = PETSC_COMM_SELF;
54340475a53SFande Kong 	}
54440475a53SFande Kong 	/*get sub-matrix data*/
5451b536640SFande Kong 	sxadj=0; sadjncy=0; svalues=0;
54685374dbcSFande Kong     ierr = MatGetSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr);
5475495542aSFande Kong     ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr);
5485495542aSFande Kong     ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr);
54985374dbcSFande Kong     ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr);
5505495542aSFande Kong     ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_n);CHKERRQ(ierr);
55185374dbcSFande Kong     ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr);
55285374dbcSFande Kong     ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr);
5535495542aSFande Kong     ierr = PetscMemcpy(indices+irow_n,icol_indices,sizeof(PetscInt)*icol_n);CHKERRQ(ierr);
55485374dbcSFande Kong     ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr);
5555495542aSFande Kong     nindx = irow_n+icol_n;
5565495542aSFande Kong     ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr);
557a1f6d831SFande Kong     /* renumber columns */
5585495542aSFande Kong     for(j=0; j<irow_n; j++){
559a1f6d831SFande Kong       for(k=sxadj[j]; k<sxadj[j+1]; k++){
5605495542aSFande Kong     	ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr);
56185374dbcSFande Kong #if PETSC_USE_DEBUG
56285374dbcSFande Kong     	if(loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %d \n",sadjncy[k]);
56385374dbcSFande Kong #endif
56485374dbcSFande Kong         sadjncy[k] = loc;
56585374dbcSFande Kong       }
56685374dbcSFande Kong     }
56785374dbcSFande Kong     if(scall==MAT_INITIAL_MATRIX){
56840475a53SFande Kong       ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr);
56985374dbcSFande Kong     }else{
57085374dbcSFande Kong        Mat                sadj = *(submat[i]);
57127157e21SSatish Balay        Mat_MPIAdj         *sa  = (Mat_MPIAdj*)((sadj)->data);
57240475a53SFande Kong        ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr);
57340475a53SFande Kong        ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRQ(ierr);
57440475a53SFande Kong        if(issame != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"submatrix  must have the same comm as the col index set\n");
5755495542aSFande Kong        ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_n+1));CHKERRQ(ierr);
5765495542aSFande Kong        ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);
5775495542aSFande Kong        if(svalues){ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);}
57885374dbcSFande Kong        ierr = PetscFree(sxadj);CHKERRQ(ierr);
57985374dbcSFande Kong        ierr = PetscFree(sadjncy);CHKERRQ(ierr);
5805495542aSFande Kong        if(svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);}
58185374dbcSFande Kong     }
58285374dbcSFande Kong   }
5835495542aSFande Kong   ierr = PetscFree(indices);CHKERRQ(ierr);
58485374dbcSFande Kong   PetscFunctionReturn(0);
58585374dbcSFande Kong }
58685374dbcSFande Kong 
58785374dbcSFande Kong 
58885374dbcSFande Kong #undef __FUNCT__
58985374dbcSFande Kong #define __FUNCT__ "MatGetSubMatrix_MPIAdj_data"
5905495542aSFande Kong /*
5915495542aSFande Kong  * The interface should be easy to use for both MatGetSubMatrix (parallel sub-matrix) and MatGetSubMatrices (sequential sub-matrices)
5925495542aSFande Kong  * */
5935495542aSFande Kong static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS irows, IS icols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values)
594841d17a1SFande Kong {
5955495542aSFande Kong   PetscInt        	 nlrows_is,icols_n,i,j,nroots,nleaves,owner,rlocalindex,*ncols_send,*ncols_recv;
5965495542aSFande Kong   PetscInt           nlrows_mat,*adjncy_recv,Ncols_recv,Ncols_send,*xadj_recv,*values_recv;
5975495542aSFande Kong   PetscInt          *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues,isvalue;
5985495542aSFande Kong   const PetscInt    *irows_indices,*icols_indices,*xadj, *adjncy;
599841d17a1SFande Kong   Mat_MPIAdj        *a = (Mat_MPIAdj*)adj->data;
600841d17a1SFande Kong   PetscLayout        rmap;
601841d17a1SFande Kong   MPI_Comm           comm;
602841d17a1SFande Kong   PetscSF            sf;
603841d17a1SFande Kong   PetscSFNode       *iremote;
604841d17a1SFande Kong   PetscBool          done;
605841d17a1SFande Kong   PetscErrorCode     ierr;
606841d17a1SFande Kong 
607841d17a1SFande Kong   PetscFunctionBegin;
608841d17a1SFande Kong   /* communicator */
609841d17a1SFande Kong   ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr);
610841d17a1SFande Kong   /* Layouts */
611841d17a1SFande Kong   ierr = MatGetLayouts(adj,&rmap,PETSC_NULL);CHKERRQ(ierr);
612841d17a1SFande Kong   /* get rows information */
6135495542aSFande Kong   ierr = ISGetLocalSize(irows,&nlrows_is);CHKERRQ(ierr);
6145495542aSFande Kong   ierr = ISGetIndices(irows,&irows_indices);CHKERRQ(ierr);
6155495542aSFande Kong   ierr = PetscCalloc1(nlrows_is,&iremote);CHKERRQ(ierr);
6165495542aSFande Kong   /* construct sf graph*/
6175495542aSFande Kong   nleaves = nlrows_is;
6185495542aSFande Kong   for(i=0; i<nlrows_is; i++){
6192253a2e4SFande Kong 	owner = -1;
6202253a2e4SFande Kong 	rlocalindex = -1;
6215495542aSFande Kong     ierr = PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex);CHKERRQ(ierr);
622841d17a1SFande Kong     iremote[i].rank  = owner;
6235495542aSFande Kong     iremote[i].index = rlocalindex;
624841d17a1SFande Kong   }
6255495542aSFande Kong   ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr);
6265495542aSFande Kong   ierr = PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv);CHKERRQ(ierr);
6275495542aSFande Kong   nroots = nlrows_mat;
6285495542aSFande Kong   for(i=0; i<nlrows_mat; i++){
629841d17a1SFande Kong 	ncols_send[i] = xadj[i+1]-xadj[i];
630841d17a1SFande Kong   }
631841d17a1SFande Kong   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
63285374dbcSFande Kong   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
633841d17a1SFande Kong   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
634841d17a1SFande Kong   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
635841d17a1SFande Kong   ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr);
636841d17a1SFande Kong   ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr);
637841d17a1SFande Kong   ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr);
638841d17a1SFande Kong   ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr);
639841d17a1SFande Kong   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
640841d17a1SFande Kong   Ncols_recv =0;
6415495542aSFande Kong   for(i=0; i<nlrows_is; i++){
642841d17a1SFande Kong 	 Ncols_recv             += ncols_recv[i];
64385374dbcSFande Kong 	 ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i];
644841d17a1SFande Kong   }
645841d17a1SFande Kong   Ncols_send = 0;
6465495542aSFande Kong   for(i=0; i<nlrows_mat; i++){
647841d17a1SFande Kong 	Ncols_send += ncols_send[i];
648841d17a1SFande Kong   }
649841d17a1SFande Kong   ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr);
6505495542aSFande Kong   ierr = PetscCalloc1(Ncols_recv,&adjncy_recv);CHKERRQ(ierr);
651841d17a1SFande Kong   nleaves = Ncols_recv;
6525495542aSFande Kong   Ncols_recv = 0;
6535495542aSFande Kong   for(i=0; i<nlrows_is; i++){
6545495542aSFande Kong     ierr = PetscLayoutFindOwner(rmap,irows_indices[i],&owner);CHKERRQ(ierr);
655841d17a1SFande Kong     for(j=0; j<ncols_recv[i]; j++){
6565495542aSFande Kong       iremote[Ncols_recv].rank    = owner;
6575495542aSFande Kong       iremote[Ncols_recv++].index = xadj_recv[i]+j;
658841d17a1SFande Kong     }
659841d17a1SFande Kong   }
6605495542aSFande Kong   ierr = ISRestoreIndices(irows,&irows_indices);CHKERRQ(ierr);
661a1f6d831SFande Kong   /*if we need to deal with edge weights ???*/
6625495542aSFande Kong   if(a->values){isvalue=1;}else{isvalue=0;}
6635495542aSFande Kong   /*involve a global communication */
664*b2566f29SBarry Smith   /*ierr = MPIU_Allreduce(&isvalue,&isvalue,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);*/
6655495542aSFande Kong   if(isvalue){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);}
666841d17a1SFande Kong   nroots = Ncols_send;
667841d17a1SFande Kong   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
66885374dbcSFande Kong   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
669841d17a1SFande Kong   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
670841d17a1SFande Kong   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
671841d17a1SFande Kong   ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr);
672841d17a1SFande Kong   ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr);
6735495542aSFande Kong   if(isvalue){
674841d17a1SFande Kong 	ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr);
675841d17a1SFande Kong 	ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr);
6765495542aSFande Kong   }
677841d17a1SFande Kong   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
6785495542aSFande Kong   ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr);
6795495542aSFande Kong   ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr);
6805495542aSFande Kong   ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr);
68185374dbcSFande Kong   rnclos = 0;
6825495542aSFande Kong   for(i=0; i<nlrows_is; i++){
68385374dbcSFande Kong     for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){
6845495542aSFande Kong       ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr);
68585374dbcSFande Kong       if(loc<0){
68685374dbcSFande Kong         adjncy_recv[j] = -1;
6875495542aSFande Kong         if(isvalue) values_recv[j] = -1;
68885374dbcSFande Kong         ncols_recv[i]--;
68985374dbcSFande Kong       }else{
69085374dbcSFande Kong     	rnclos++;
69185374dbcSFande Kong       }
69285374dbcSFande Kong     }
69385374dbcSFande Kong   }
6945495542aSFande Kong   ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr);
69585374dbcSFande Kong   ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr);
6965495542aSFande Kong   if(isvalue) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);}
6975495542aSFande Kong   ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr);
69885374dbcSFande Kong   rnclos = 0;
6995495542aSFande Kong   for(i=0; i<nlrows_is; i++){
70085374dbcSFande Kong 	for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){
70185374dbcSFande Kong 	  if(adjncy_recv[j]<0) continue;
70285374dbcSFande Kong 	  sadjncy[rnclos] = adjncy_recv[j];
7035495542aSFande Kong 	  if(isvalue) svalues[rnclos] = values_recv[j];
70485374dbcSFande Kong 	  rnclos++;
70585374dbcSFande Kong 	}
70685374dbcSFande Kong   }
7075495542aSFande Kong   for(i=0; i<nlrows_is; i++){
70885374dbcSFande Kong 	sxadj[i+1] = sxadj[i]+ncols_recv[i];
70985374dbcSFande Kong   }
71085374dbcSFande Kong   if(sadj_xadj)  { *sadj_xadj = sxadj;}else    { ierr = PetscFree(sxadj);CHKERRQ(ierr);}
71185374dbcSFande Kong   if(sadj_adjncy){ *sadj_adjncy = sadjncy;}else{ ierr = PetscFree(sadjncy);CHKERRQ(ierr);}
7125495542aSFande Kong   if(sadj_values){
7135495542aSFande Kong 	if(isvalue) *sadj_values = svalues; else *sadj_values=0;
7145495542aSFande Kong   }else{
7155495542aSFande Kong 	if(isvalue) {ierr = PetscFree(svalues);CHKERRQ(ierr);}
7165495542aSFande Kong   }
71785374dbcSFande Kong   ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr);
7185495542aSFande Kong   ierr = PetscFree(adjncy_recv);CHKERRQ(ierr);
7195495542aSFande Kong   if(isvalue) {ierr = PetscFree(values_recv);CHKERRQ(ierr);}
720841d17a1SFande Kong   PetscFunctionReturn(0);
721841d17a1SFande Kong }
722841d17a1SFande Kong 
723841d17a1SFande Kong 
724841d17a1SFande Kong #undef __FUNCT__
7259a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj"
726f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
7279a3665c6SJed Brown {
7289a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
7299a3665c6SJed Brown   PetscErrorCode ierr;
7309a3665c6SJed Brown   const PetscInt *ranges;
7319a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
7329a3665c6SJed Brown   MPI_Group      agroup,bgroup;
7339a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
7349a3665c6SJed Brown 
7359a3665c6SJed Brown   PetscFunctionBegin;
7360298fd71SBarry Smith   *B    = NULL;
737ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
7389a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
7399a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
7409a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
7419a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
7429a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
7439a3665c6SJed Brown   }
7449a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
7459a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
7469a3665c6SJed Brown     *B   = A;
7479a3665c6SJed Brown     PetscFunctionReturn(0);
7489a3665c6SJed Brown   }
7499a3665c6SJed Brown 
750785e854fSJed Brown   ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr);
7519a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
7529a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
7539a3665c6SJed Brown   }
7549a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
7559a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
7569a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
7579a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
7589a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
7599a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
7609a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
7619a3665c6SJed Brown     PetscInt   m,N;
7629a3665c6SJed Brown     Mat_MPIAdj *b;
7630298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
7640298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
7659a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
7669a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
7679a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
7689a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
7699a3665c6SJed Brown   }
7709a3665c6SJed Brown   PetscFunctionReturn(0);
7719a3665c6SJed Brown }
7729a3665c6SJed Brown 
7739a3665c6SJed Brown #undef __FUNCT__
7749a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat"
7759a3665c6SJed Brown /*@
7769a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
7779a3665c6SJed Brown 
7789a3665c6SJed Brown    Collective
7799a3665c6SJed Brown 
7809a3665c6SJed Brown    Input Arguments:
7819a3665c6SJed Brown .  A - original MPIAdj matrix
7829a3665c6SJed Brown 
7839a3665c6SJed Brown    Output Arguments:
7840298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
7859a3665c6SJed Brown 
7869a3665c6SJed Brown    Level: developer
7879a3665c6SJed Brown 
7889a3665c6SJed Brown    Note:
7899a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
7909a3665c6SJed Brown 
7919a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
7929a3665c6SJed Brown 
7939a3665c6SJed Brown .seealso: MatCreateMPIAdj()
7949a3665c6SJed Brown @*/
7959a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
7969a3665c6SJed Brown {
7979a3665c6SJed Brown   PetscErrorCode ierr;
7989a3665c6SJed Brown 
7999a3665c6SJed Brown   PetscFunctionBegin;
8009a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8019a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
8029a3665c6SJed Brown   PetscFunctionReturn(0);
8039a3665c6SJed Brown }
8049a3665c6SJed Brown 
8050bad9183SKris Buschelman /*MC
806fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
8070bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
8080bad9183SKris Buschelman 
8090bad9183SKris Buschelman   Level: beginner
8100bad9183SKris Buschelman 
8110bad9183SKris Buschelman .seealso: MatCreateMPIAdj
8120bad9183SKris Buschelman M*/
8130bad9183SKris Buschelman 
8144a2ae208SSatish Balay #undef __FUNCT__
8154a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
8168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
817273d9f13SBarry Smith {
818273d9f13SBarry Smith   Mat_MPIAdj     *b;
8196849ba73SBarry Smith   PetscErrorCode ierr;
820273d9f13SBarry Smith 
821273d9f13SBarry Smith   PetscFunctionBegin;
822b00a9115SJed Brown   ierr         = PetscNewLog(B,&b);CHKERRQ(ierr);
823b0a32e0cSBarry Smith   B->data      = (void*)b;
824273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
825273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
826273d9f13SBarry Smith 
827bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
828bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
82917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
830273d9f13SBarry Smith   PetscFunctionReturn(0);
831273d9f13SBarry Smith }
832273d9f13SBarry Smith 
8334a2ae208SSatish Balay #undef __FUNCT__
8344a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
835a23d5eceSKris Buschelman /*@C
836a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
837a23d5eceSKris Buschelman 
8383f9fe445SBarry Smith    Logically Collective on MPI_Comm
839a23d5eceSKris Buschelman 
840a23d5eceSKris Buschelman    Input Parameters:
841a23d5eceSKris Buschelman +  A - the matrix
842a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
843a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
844a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
845a23d5eceSKris Buschelman -  values - [optional] edge weights
846a23d5eceSKris Buschelman 
847a23d5eceSKris Buschelman    Level: intermediate
848a23d5eceSKris Buschelman 
849a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
850a23d5eceSKris Buschelman @*/
8517087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
852273d9f13SBarry Smith {
8534ac538c5SBarry Smith   PetscErrorCode ierr;
854273d9f13SBarry Smith 
855273d9f13SBarry Smith   PetscFunctionBegin;
8564ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
857273d9f13SBarry Smith   PetscFunctionReturn(0);
858273d9f13SBarry Smith }
859273d9f13SBarry Smith 
8604a2ae208SSatish Balay #undef __FUNCT__
8614a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
862b97cf49bSBarry Smith /*@C
8633eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
864b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
865b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
866b97cf49bSBarry Smith 
867ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
868ef5ee4d1SLois Curfman McInnes 
869b97cf49bSBarry Smith    Input Parameters:
870c2e958feSBarry Smith +  comm - MPI communicator
8710752156aSBarry Smith .  m - number of local rows
87258ec18a5SLisandro Dalcin .  N - number of global columns
873b97cf49bSBarry Smith .  i - the indices into j for the start of each row
874329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
875ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
876329f5518SBarry Smith -  values -[optional] edge weights
877b97cf49bSBarry Smith 
878b97cf49bSBarry Smith    Output Parameter:
879b97cf49bSBarry Smith .  A - the matrix
880b97cf49bSBarry Smith 
88115091d37SBarry Smith    Level: intermediate
88215091d37SBarry Smith 
8834bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
8844bc6d8c0SBarry Smith    MatSetValues().
885c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
886ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
8871198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
888327e1a73SBarry Smith    Should not include the matrix diagonals.
889b97cf49bSBarry Smith 
890e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
891ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
892ededeb1aSvictorle 
893ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
894b97cf49bSBarry Smith 
895e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
896b97cf49bSBarry Smith @*/
8977087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
898b97cf49bSBarry Smith {
899dfbe8321SBarry Smith   PetscErrorCode ierr;
900b97cf49bSBarry Smith 
901433994e6SBarry Smith   PetscFunctionBegin;
902f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
90358ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
904273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
905273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
9063a40ed3dSBarry Smith   PetscFunctionReturn(0);
907b97cf49bSBarry Smith }
908b97cf49bSBarry Smith 
909c2e958feSBarry Smith 
910433994e6SBarry Smith 
911433994e6SBarry Smith 
912433994e6SBarry Smith 
913433994e6SBarry Smith 
914433994e6SBarry Smith 
915