xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision f9426fe092dba0ba2fdf65dfec8d938c4b10a31c)
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*/
6b97cf49bSBarry Smith 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
9dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
10b97cf49bSBarry Smith {
113eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
12dfbe8321SBarry Smith   PetscErrorCode    ierr;
13d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n;
142dcb1b2aSMatthew Knepley   const char        *name;
15ce1411ecSBarry Smith   PetscViewerFormat format;
16b97cf49bSBarry Smith 
17433994e6SBarry Smith   PetscFunctionBegin;
183a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
19b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
20fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
213a40ed3dSBarry Smith     PetscFunctionReturn(0);
22fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
23ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported");
240752156aSBarry Smith   } else {
25d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
267b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);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);
367b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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   }
78bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
79dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
80bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
81bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
823a40ed3dSBarry Smith   PetscFunctionReturn(0);
83b97cf49bSBarry Smith }
84b97cf49bSBarry Smith 
854a2ae208SSatish Balay #undef __FUNCT__
864a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
87ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
88b97cf49bSBarry Smith {
893eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
9063ba0a88SBarry Smith   PetscErrorCode ierr;
91b97cf49bSBarry Smith 
92433994e6SBarry Smith   PetscFunctionBegin;
9312c028f9SKris Buschelman   switch (op) {
9477e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9512c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
969a4540c5SBarry Smith   case MAT_HERMITIAN:
974e0d8c25SBarry Smith     a->symmetric = flg;
989a4540c5SBarry Smith     break;
999a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1009a4540c5SBarry Smith     break;
10112c028f9SKris Buschelman   default:
102290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
10312c028f9SKris Buschelman     break;
104b97cf49bSBarry Smith   }
1053a40ed3dSBarry Smith   PetscFunctionReturn(0);
106b97cf49bSBarry Smith }
107b97cf49bSBarry Smith 
108b97cf49bSBarry Smith 
109b97cf49bSBarry Smith /*
110b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
111b97cf49bSBarry Smith */
112b97cf49bSBarry Smith 
1134a2ae208SSatish Balay #undef __FUNCT__
1144a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
115dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
116b97cf49bSBarry Smith {
1173eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1186849ba73SBarry Smith   PetscErrorCode ierr;
119d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n;
120b97cf49bSBarry Smith 
121433994e6SBarry Smith   PetscFunctionBegin;
12209f38230SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
1233bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,m*sizeof(PetscInt));CHKERRQ(ierr);
124d0f46423SBarry Smith   for (i=0; i<A->rmap->n; i++) {
125b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
126b97cf49bSBarry Smith       if (a->j[j] == i) {
12709f38230SBarry Smith         a->diag[i] = j;
128b97cf49bSBarry Smith         break;
129b97cf49bSBarry Smith       }
130b97cf49bSBarry Smith     }
131b97cf49bSBarry Smith   }
1323a40ed3dSBarry Smith   PetscFunctionReturn(0);
133b97cf49bSBarry Smith }
134b97cf49bSBarry Smith 
1354a2ae208SSatish Balay #undef __FUNCT__
1364a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
137b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
138b97cf49bSBarry Smith {
1393eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
140b24ad042SBarry Smith   PetscInt   *itmp;
141b97cf49bSBarry Smith 
142433994e6SBarry Smith   PetscFunctionBegin;
143d0f46423SBarry Smith   row -= A->rmap->rstart;
1440752156aSBarry Smith 
145e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
146b97cf49bSBarry Smith 
147b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
1480298fd71SBarry Smith   if (v) *v = NULL;
149b97cf49bSBarry Smith   if (idx) {
150b97cf49bSBarry Smith     itmp = a->j + a->i[row];
151b97cf49bSBarry Smith     if (*nz) {
152b97cf49bSBarry Smith       *idx = itmp;
153b97cf49bSBarry Smith     }
154b97cf49bSBarry Smith     else *idx = 0;
155b97cf49bSBarry Smith   }
1563a40ed3dSBarry Smith   PetscFunctionReturn(0);
157b97cf49bSBarry Smith }
158b97cf49bSBarry Smith 
1594a2ae208SSatish Balay #undef __FUNCT__
1604a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
161b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
162b97cf49bSBarry Smith {
163433994e6SBarry Smith   PetscFunctionBegin;
1643a40ed3dSBarry Smith   PetscFunctionReturn(0);
165b97cf49bSBarry Smith }
166b97cf49bSBarry Smith 
1674a2ae208SSatish Balay #undef __FUNCT__
1684a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
169ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
170b97cf49bSBarry Smith {
1713eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
172dfbe8321SBarry Smith   PetscErrorCode ierr;
173ace3abfcSBarry Smith   PetscBool      flag;
174b97cf49bSBarry Smith 
175433994e6SBarry Smith   PetscFunctionBegin;
176b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
177d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1780f5bd95cSBarry Smith     flag = PETSC_FALSE;
179b97cf49bSBarry Smith   }
180b97cf49bSBarry Smith 
181b97cf49bSBarry Smith   /* if the a->i are the same */
182d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
183b97cf49bSBarry Smith 
184b97cf49bSBarry Smith   /* if a->j are the same */
185b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
186b97cf49bSBarry Smith 
187ce94432eSBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
189b97cf49bSBarry Smith }
190b97cf49bSBarry Smith 
1914a2ae208SSatish Balay #undef __FUNCT__
1924a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
1931a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
1946cd91f64SBarry Smith {
195b24ad042SBarry Smith   PetscInt   i;
1966cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
1971a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
1986cd91f64SBarry Smith 
1996cd91f64SBarry Smith   PetscFunctionBegin;
200d0f46423SBarry Smith   *m    = A->rmap->n;
2016cd91f64SBarry Smith   *ia   = a->i;
2026cd91f64SBarry Smith   *ja   = a->j;
2036cd91f64SBarry Smith   *done = PETSC_TRUE;
204d42735a1SBarry Smith   if (oshift) {
205d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
206d42735a1SBarry Smith       (*ja)[i]++;
207d42735a1SBarry Smith     }
208d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
209d42735a1SBarry Smith   }
210d42735a1SBarry Smith   PetscFunctionReturn(0);
211d42735a1SBarry Smith }
212d42735a1SBarry Smith 
2134a2ae208SSatish Balay #undef __FUNCT__
2144a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
2151a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
216d42735a1SBarry Smith {
217b24ad042SBarry Smith   PetscInt   i;
218d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2191a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
220d42735a1SBarry Smith 
221d42735a1SBarry Smith   PetscFunctionBegin;
222e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
223e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
224d42735a1SBarry Smith   if (oshift) {
225d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
226d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
227d42735a1SBarry Smith       (*ja)[i]--;
228d42735a1SBarry Smith     }
229d42735a1SBarry Smith   }
2306cd91f64SBarry Smith   PetscFunctionReturn(0);
2316cd91f64SBarry Smith }
232b97cf49bSBarry Smith 
23317667f90SBarry Smith #undef __FUNCT__
23417667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
23519fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
23617667f90SBarry Smith {
23717667f90SBarry Smith   Mat               B;
23817667f90SBarry Smith   PetscErrorCode    ierr;
23917667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
24017667f90SBarry Smith   const PetscInt    *rj;
24117667f90SBarry Smith   const PetscScalar *ra;
24217667f90SBarry Smith   MPI_Comm          comm;
24317667f90SBarry Smith 
24417667f90SBarry Smith   PetscFunctionBegin;
2450298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
2460298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
2470298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
24817667f90SBarry Smith 
24917667f90SBarry Smith   /* count the number of nonzeros per row */
25017667f90SBarry Smith   for (i=0; i<m; i++) {
2510298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
2525ee9ba1cSJed Brown     for (j=0; j<len; j++) {
2535ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
2545ee9ba1cSJed Brown     }
25517667f90SBarry Smith     nzeros += len;
2560f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
25717667f90SBarry Smith   }
25817667f90SBarry Smith 
25917667f90SBarry Smith   /* malloc space for nonzeros */
26017667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr);
26117667f90SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
26217667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr);
26317667f90SBarry Smith 
26417667f90SBarry Smith   nzeros = 0;
26517667f90SBarry Smith   ia[0]  = 0;
26617667f90SBarry Smith   for (i=0; i<m; i++) {
26717667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
26817667f90SBarry Smith     cnt  = 0;
26917667f90SBarry Smith     for (j=0; j<len; j++) {
2705ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
27117667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27217667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
27317667f90SBarry Smith       }
2745ee9ba1cSJed Brown     }
27517667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27617667f90SBarry Smith     nzeros += cnt;
27717667f90SBarry Smith     ia[i+1] = nzeros;
27817667f90SBarry Smith   }
27917667f90SBarry Smith 
28017667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
28117667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
28258ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28317667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28417667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
28517667f90SBarry Smith 
28617667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
28717667f90SBarry Smith     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
28817667f90SBarry Smith   } else {
28917667f90SBarry Smith     *newmat = B;
29017667f90SBarry Smith   }
29117667f90SBarry Smith   PetscFunctionReturn(0);
29217667f90SBarry Smith }
29317667f90SBarry Smith 
294b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2950b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2963eda8832SBarry Smith                                        MatGetRow_MPIAdj,
2973eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
298b97cf49bSBarry Smith                                        0,
29997304618SKris Buschelman                                 /* 4*/ 0,
300b97cf49bSBarry Smith                                        0,
301b97cf49bSBarry Smith                                        0,
302b97cf49bSBarry Smith                                        0,
3030b3b1487SBarry Smith                                        0,
3040b3b1487SBarry Smith                                        0,
30597304618SKris Buschelman                                 /*10*/ 0,
3060b3b1487SBarry Smith                                        0,
3070b3b1487SBarry Smith                                        0,
3080b3b1487SBarry Smith                                        0,
3090b3b1487SBarry Smith                                        0,
31097304618SKris Buschelman                                 /*15*/ 0,
3113eda8832SBarry Smith                                        MatEqual_MPIAdj,
3120b3b1487SBarry Smith                                        0,
3130b3b1487SBarry Smith                                        0,
3140b3b1487SBarry Smith                                        0,
31597304618SKris Buschelman                                 /*20*/ 0,
3160b3b1487SBarry Smith                                        0,
3173eda8832SBarry Smith                                        MatSetOption_MPIAdj,
3180b3b1487SBarry Smith                                        0,
319d519adbfSMatthew Knepley                                 /*24*/ 0,
3200b3b1487SBarry Smith                                        0,
3210b3b1487SBarry Smith                                        0,
3220b3b1487SBarry Smith                                        0,
3230b3b1487SBarry Smith                                        0,
324d519adbfSMatthew Knepley                                 /*29*/ 0,
3250b3b1487SBarry Smith                                        0,
326273d9f13SBarry Smith                                        0,
327273d9f13SBarry Smith                                        0,
3280b3b1487SBarry Smith                                        0,
329d519adbfSMatthew Knepley                                 /*34*/ 0,
3300b3b1487SBarry Smith                                        0,
3310b3b1487SBarry Smith                                        0,
3320b3b1487SBarry Smith                                        0,
3330b3b1487SBarry Smith                                        0,
334d519adbfSMatthew Knepley                                 /*39*/ 0,
3350b3b1487SBarry Smith                                        0,
3360b3b1487SBarry Smith                                        0,
3370b3b1487SBarry Smith                                        0,
3380b3b1487SBarry Smith                                        0,
339d519adbfSMatthew Knepley                                 /*44*/ 0,
3400b3b1487SBarry Smith                                        0,
3410b3b1487SBarry Smith                                        0,
3420b3b1487SBarry Smith                                        0,
3430b3b1487SBarry Smith                                        0,
344d519adbfSMatthew Knepley                                 /*49*/ 0,
3456cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
346d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
347b97cf49bSBarry Smith                                        0,
348b97cf49bSBarry Smith                                        0,
349d519adbfSMatthew Knepley                                 /*54*/ 0,
350b97cf49bSBarry Smith                                        0,
3510752156aSBarry Smith                                        0,
3520752156aSBarry Smith                                        0,
3530b3b1487SBarry Smith                                        0,
354d519adbfSMatthew Knepley                                 /*59*/ 0,
355b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
356b9b97703SBarry Smith                                        MatView_MPIAdj,
35717667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
35897304618SKris Buschelman                                        0,
359d519adbfSMatthew Knepley                                 /*64*/ 0,
36097304618SKris Buschelman                                        0,
36197304618SKris Buschelman                                        0,
36297304618SKris Buschelman                                        0,
36397304618SKris Buschelman                                        0,
364d519adbfSMatthew Knepley                                 /*69*/ 0,
36597304618SKris Buschelman                                        0,
36697304618SKris Buschelman                                        0,
36797304618SKris Buschelman                                        0,
36897304618SKris Buschelman                                        0,
369d519adbfSMatthew Knepley                                 /*74*/ 0,
37097304618SKris Buschelman                                        0,
37197304618SKris Buschelman                                        0,
37297304618SKris Buschelman                                        0,
37397304618SKris Buschelman                                        0,
374d519adbfSMatthew Knepley                                 /*79*/ 0,
37597304618SKris Buschelman                                        0,
37697304618SKris Buschelman                                        0,
37797304618SKris Buschelman                                        0,
378865e5f61SKris Buschelman                                        0,
379d519adbfSMatthew Knepley                                 /*84*/ 0,
380865e5f61SKris Buschelman                                        0,
381865e5f61SKris Buschelman                                        0,
382865e5f61SKris Buschelman                                        0,
383865e5f61SKris Buschelman                                        0,
384d519adbfSMatthew Knepley                                 /*89*/ 0,
385865e5f61SKris Buschelman                                        0,
386865e5f61SKris Buschelman                                        0,
387865e5f61SKris Buschelman                                        0,
388865e5f61SKris Buschelman                                        0,
389d519adbfSMatthew Knepley                                 /*94*/ 0,
390865e5f61SKris Buschelman                                        0,
391865e5f61SKris Buschelman                                        0,
3923964eb88SJed Brown                                        0,
3933964eb88SJed Brown                                        0,
3943964eb88SJed Brown                                 /*99*/ 0,
3953964eb88SJed Brown                                        0,
3963964eb88SJed Brown                                        0,
3973964eb88SJed Brown                                        0,
3983964eb88SJed Brown                                        0,
3993964eb88SJed Brown                                /*104*/ 0,
4003964eb88SJed Brown                                        0,
4013964eb88SJed Brown                                        0,
4023964eb88SJed Brown                                        0,
4033964eb88SJed Brown                                        0,
4043964eb88SJed Brown                                /*109*/ 0,
4053964eb88SJed Brown                                        0,
4063964eb88SJed Brown                                        0,
4073964eb88SJed Brown                                        0,
4083964eb88SJed Brown                                        0,
4093964eb88SJed Brown                                /*114*/ 0,
4103964eb88SJed Brown                                        0,
4113964eb88SJed Brown                                        0,
4123964eb88SJed Brown                                        0,
4133964eb88SJed Brown                                        0,
4143964eb88SJed Brown                                /*119*/ 0,
4153964eb88SJed Brown                                        0,
4163964eb88SJed Brown                                        0,
4173964eb88SJed Brown                                        0,
4183964eb88SJed Brown                                        0,
4193964eb88SJed Brown                                /*124*/ 0,
4203964eb88SJed Brown                                        0,
4213964eb88SJed Brown                                        0,
4223964eb88SJed Brown                                        0,
4233964eb88SJed Brown                                        0,
4243964eb88SJed Brown                                /*129*/ 0,
4253964eb88SJed Brown                                        0,
4263964eb88SJed Brown                                        0,
4273964eb88SJed Brown                                        0,
4283964eb88SJed Brown                                        0,
4293964eb88SJed Brown                                /*134*/ 0,
4303964eb88SJed Brown                                        0,
4313964eb88SJed Brown                                        0,
4323964eb88SJed Brown                                        0,
4333964eb88SJed Brown                                        0,
4343964eb88SJed Brown                                /*139*/ 0,
435*f9426fe0SMark Adams                                        0,
4363964eb88SJed Brown                                        0
4373964eb88SJed Brown };
438b97cf49bSBarry Smith 
439a23d5eceSKris Buschelman #undef __FUNCT__
440a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
441f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
442a23d5eceSKris Buschelman {
443a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
444dfbe8321SBarry Smith   PetscErrorCode ierr;
4452515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
446b24ad042SBarry Smith   PetscInt ii;
447a23d5eceSKris Buschelman #endif
448a23d5eceSKris Buschelman 
449a23d5eceSKris Buschelman   PetscFunctionBegin;
45026283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
45126283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
45258ec18a5SLisandro Dalcin 
4532515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
454e32f2f54SBarry 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]);
455d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
456e02043d6SBarry 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]);
457a23d5eceSKris Buschelman   }
458d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
459e02043d6SBarry 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]);
460a23d5eceSKris Buschelman   }
461a23d5eceSKris Buschelman #endif
46258ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
463a23d5eceSKris Buschelman 
464a23d5eceSKris Buschelman   b->j      = j;
465a23d5eceSKris Buschelman   b->i      = i;
466a23d5eceSKris Buschelman   b->values = values;
467a23d5eceSKris Buschelman 
468d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
469a23d5eceSKris Buschelman   b->diag      = 0;
470a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
471a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
472a23d5eceSKris Buschelman 
473a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
474a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
475a23d5eceSKris Buschelman   PetscFunctionReturn(0);
476a23d5eceSKris Buschelman }
477b97cf49bSBarry Smith 
4789a3665c6SJed Brown #undef __FUNCT__
4799a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj"
480f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
4819a3665c6SJed Brown {
4829a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
4839a3665c6SJed Brown   PetscErrorCode ierr;
4849a3665c6SJed Brown   const PetscInt *ranges;
4859a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
4869a3665c6SJed Brown   MPI_Group      agroup,bgroup;
4879a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
4889a3665c6SJed Brown 
4899a3665c6SJed Brown   PetscFunctionBegin;
4900298fd71SBarry Smith   *B    = NULL;
491ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
4929a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
4939a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
4949a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
4959a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
4969a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
4979a3665c6SJed Brown   }
4989a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
4999a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
5009a3665c6SJed Brown     *B   = A;
5019a3665c6SJed Brown     PetscFunctionReturn(0);
5029a3665c6SJed Brown   }
5039a3665c6SJed Brown 
5049a3665c6SJed Brown   ierr = PetscMalloc(nranks*sizeof(PetscMPIInt),&ranks);CHKERRQ(ierr);
5059a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
5069a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
5079a3665c6SJed Brown   }
5089a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
5099a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
5109a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
5119a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
5129a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
5139a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
5149a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
5159a3665c6SJed Brown     PetscInt   m,N;
5169a3665c6SJed Brown     Mat_MPIAdj *b;
5170298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
5180298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
5199a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
5209a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
5219a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
5229a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
5239a3665c6SJed Brown   }
5249a3665c6SJed Brown   PetscFunctionReturn(0);
5259a3665c6SJed Brown }
5269a3665c6SJed Brown 
5279a3665c6SJed Brown #undef __FUNCT__
5289a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat"
5299a3665c6SJed Brown /*@
5309a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
5319a3665c6SJed Brown 
5329a3665c6SJed Brown    Collective
5339a3665c6SJed Brown 
5349a3665c6SJed Brown    Input Arguments:
5359a3665c6SJed Brown .  A - original MPIAdj matrix
5369a3665c6SJed Brown 
5379a3665c6SJed Brown    Output Arguments:
5380298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
5399a3665c6SJed Brown 
5409a3665c6SJed Brown    Level: developer
5419a3665c6SJed Brown 
5429a3665c6SJed Brown    Note:
5439a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
5449a3665c6SJed Brown 
5459a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
5469a3665c6SJed Brown 
5479a3665c6SJed Brown .seealso: MatCreateMPIAdj()
5489a3665c6SJed Brown @*/
5499a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
5509a3665c6SJed Brown {
5519a3665c6SJed Brown   PetscErrorCode ierr;
5529a3665c6SJed Brown 
5539a3665c6SJed Brown   PetscFunctionBegin;
5549a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
5559a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
5569a3665c6SJed Brown   PetscFunctionReturn(0);
5579a3665c6SJed Brown }
5589a3665c6SJed Brown 
5590bad9183SKris Buschelman /*MC
560fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
5610bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
5620bad9183SKris Buschelman 
5630bad9183SKris Buschelman   Level: beginner
5640bad9183SKris Buschelman 
5650bad9183SKris Buschelman .seealso: MatCreateMPIAdj
5660bad9183SKris Buschelman M*/
5670bad9183SKris Buschelman 
5684a2ae208SSatish Balay #undef __FUNCT__
5694a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
5708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
571273d9f13SBarry Smith {
572273d9f13SBarry Smith   Mat_MPIAdj     *b;
5736849ba73SBarry Smith   PetscErrorCode ierr;
574273d9f13SBarry Smith 
575273d9f13SBarry Smith   PetscFunctionBegin;
57638f2d2fdSLisandro Dalcin   ierr         = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr);
577b0a32e0cSBarry Smith   B->data      = (void*)b;
578273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
579273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
580273d9f13SBarry Smith 
581bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
582bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
58317667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
584273d9f13SBarry Smith   PetscFunctionReturn(0);
585273d9f13SBarry Smith }
586273d9f13SBarry Smith 
5874a2ae208SSatish Balay #undef __FUNCT__
5884a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
589a23d5eceSKris Buschelman /*@C
590a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
591a23d5eceSKris Buschelman 
5923f9fe445SBarry Smith    Logically Collective on MPI_Comm
593a23d5eceSKris Buschelman 
594a23d5eceSKris Buschelman    Input Parameters:
595a23d5eceSKris Buschelman +  A - the matrix
596a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
597a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
598a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
599a23d5eceSKris Buschelman -  values - [optional] edge weights
600a23d5eceSKris Buschelman 
601a23d5eceSKris Buschelman    Level: intermediate
602a23d5eceSKris Buschelman 
603a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
604a23d5eceSKris Buschelman @*/
6057087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
606273d9f13SBarry Smith {
6074ac538c5SBarry Smith   PetscErrorCode ierr;
608273d9f13SBarry Smith 
609273d9f13SBarry Smith   PetscFunctionBegin;
6104ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
611273d9f13SBarry Smith   PetscFunctionReturn(0);
612273d9f13SBarry Smith }
613273d9f13SBarry Smith 
6144a2ae208SSatish Balay #undef __FUNCT__
6154a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
616b97cf49bSBarry Smith /*@C
6173eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
618b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
619b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
620b97cf49bSBarry Smith 
621ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
622ef5ee4d1SLois Curfman McInnes 
623b97cf49bSBarry Smith    Input Parameters:
624c2e958feSBarry Smith +  comm - MPI communicator
6250752156aSBarry Smith .  m - number of local rows
62658ec18a5SLisandro Dalcin .  N - number of global columns
627b97cf49bSBarry Smith .  i - the indices into j for the start of each row
628329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
629ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
630329f5518SBarry Smith -  values -[optional] edge weights
631b97cf49bSBarry Smith 
632b97cf49bSBarry Smith    Output Parameter:
633b97cf49bSBarry Smith .  A - the matrix
634b97cf49bSBarry Smith 
63515091d37SBarry Smith    Level: intermediate
63615091d37SBarry Smith 
6374bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
6384bc6d8c0SBarry Smith    MatSetValues().
639c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
640ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
6411198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
642327e1a73SBarry Smith    Should not include the matrix diagonals.
643b97cf49bSBarry Smith 
644e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
645ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
646ededeb1aSvictorle 
647ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
648b97cf49bSBarry Smith 
649e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
650b97cf49bSBarry Smith @*/
6517087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
652b97cf49bSBarry Smith {
653dfbe8321SBarry Smith   PetscErrorCode ierr;
654b97cf49bSBarry Smith 
655433994e6SBarry Smith   PetscFunctionBegin;
656f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
65758ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
658273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
659273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
661b97cf49bSBarry Smith }
662b97cf49bSBarry Smith 
663c2e958feSBarry Smith 
664433994e6SBarry Smith 
665433994e6SBarry Smith 
666433994e6SBarry Smith 
667433994e6SBarry Smith 
668433994e6SBarry Smith 
669