xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision e3c5b3bad8a8c55ce6fdc25c1dde664f52a31a1c)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3b97cf49bSBarry Smith /*
4b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
5b97cf49bSBarry Smith */
67c4f633dSBarry Smith #include "../src/mat/impls/adj/mpi/mpiadj.h"
7b97cf49bSBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
10dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
11b97cf49bSBarry Smith {
123eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
13dfbe8321SBarry Smith   PetscErrorCode    ierr;
14d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n;
152dcb1b2aSMatthew Knepley   const char        *name;
16ce1411ecSBarry Smith   PetscViewerFormat format;
17b97cf49bSBarry Smith 
18433994e6SBarry Smith   PetscFunctionBegin;
193a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
20b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
21fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
223a40ed3dSBarry Smith     PetscFunctionReturn(0);
23fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
24*e3c5b3baSBarry Smith     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MATLAB format not supported");
250752156aSBarry Smith   } else {
26d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);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);
35b97cf49bSBarry Smith   }
36b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
373a40ed3dSBarry Smith   PetscFunctionReturn(0);
38b97cf49bSBarry Smith }
39b97cf49bSBarry Smith 
404a2ae208SSatish Balay #undef __FUNCT__
414a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
42dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
43b97cf49bSBarry Smith {
44dfbe8321SBarry Smith   PetscErrorCode ierr;
45ace3abfcSBarry Smith   PetscBool      iascii;
46b97cf49bSBarry Smith 
47433994e6SBarry Smith   PetscFunctionBegin;
482692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4932077d6dSBarry Smith   if (iascii) {
503eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
515cd90555SBarry Smith   } else {
52e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
53b97cf49bSBarry Smith   }
543a40ed3dSBarry Smith   PetscFunctionReturn(0);
55b97cf49bSBarry Smith }
56b97cf49bSBarry Smith 
574a2ae208SSatish Balay #undef __FUNCT__
584a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
59dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat)
60b97cf49bSBarry Smith {
613eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
62dfbe8321SBarry Smith   PetscErrorCode ierr;
6394d884c6SBarry Smith 
6494d884c6SBarry Smith   PetscFunctionBegin;
65aa482453SBarry Smith #if defined(PETSC_USE_LOG)
66d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
67b97cf49bSBarry Smith #endif
6805b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
690400557aSBarry Smith   if (a->freeaij) {
7014196391SBarry Smith     if (a->freeaijwithfree) {
7114196391SBarry Smith       if (a->i) free(a->i);
7214196391SBarry Smith       if (a->j) free(a->j);
7314196391SBarry Smith     } else {
74606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
75606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
7605b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
7714196391SBarry Smith     }
780400557aSBarry Smith   }
791198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
80dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
81901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_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);
12309f38230SBarry Smith   ierr = PetscLogObjectMemory(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];
148b97cf49bSBarry Smith   if (v) *v = PETSC_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 
1877adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
189b97cf49bSBarry Smith }
190b97cf49bSBarry Smith 
1914a2ae208SSatish Balay #undef __FUNCT__
1924a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
193ace3abfcSBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool  *done)
1946cd91f64SBarry Smith {
195b24ad042SBarry Smith   PetscInt       i;
1966cd91f64SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
1976cd91f64SBarry Smith 
1986cd91f64SBarry Smith   PetscFunctionBegin;
199d0f46423SBarry Smith   *m    = A->rmap->n;
2006cd91f64SBarry Smith   *ia   = a->i;
2016cd91f64SBarry Smith   *ja   = a->j;
2026cd91f64SBarry Smith   *done = PETSC_TRUE;
203d42735a1SBarry Smith   if (oshift) {
204d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
205d42735a1SBarry Smith       (*ja)[i]++;
206d42735a1SBarry Smith     }
207d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
208d42735a1SBarry Smith   }
209d42735a1SBarry Smith   PetscFunctionReturn(0);
210d42735a1SBarry Smith }
211d42735a1SBarry Smith 
2124a2ae208SSatish Balay #undef __FUNCT__
2134a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
214ace3abfcSBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool  *done)
215d42735a1SBarry Smith {
216b24ad042SBarry Smith   PetscInt   i;
217d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
218d42735a1SBarry Smith 
219d42735a1SBarry Smith   PetscFunctionBegin;
220e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
221e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
222d42735a1SBarry Smith   if (oshift) {
223d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
224d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
225d42735a1SBarry Smith       (*ja)[i]--;
226d42735a1SBarry Smith     }
227d42735a1SBarry Smith   }
2286cd91f64SBarry Smith   PetscFunctionReturn(0);
2296cd91f64SBarry Smith }
230b97cf49bSBarry Smith 
23117667f90SBarry Smith #undef __FUNCT__
23217667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
2337087cfbeSBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
23417667f90SBarry Smith {
23517667f90SBarry Smith   Mat               B;
23617667f90SBarry Smith   PetscErrorCode    ierr;
23717667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
23817667f90SBarry Smith   const PetscInt    *rj;
23917667f90SBarry Smith   const PetscScalar *ra;
24017667f90SBarry Smith   MPI_Comm          comm;
24117667f90SBarry Smith 
24217667f90SBarry Smith   PetscFunctionBegin;
24317667f90SBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
24417667f90SBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
24517667f90SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
24617667f90SBarry Smith 
24717667f90SBarry Smith   /* count the number of nonzeros per row */
24817667f90SBarry Smith   for (i=0; i<m; i++) {
24917667f90SBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
25017667f90SBarry Smith     for (j=0; j<len; j++) {
25117667f90SBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
25217667f90SBarry Smith     }
25317667f90SBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
25417667f90SBarry Smith     nzeros += len;
25517667f90SBarry Smith   }
25617667f90SBarry Smith 
25717667f90SBarry Smith   /* malloc space for nonzeros */
25817667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr);
25917667f90SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
26017667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr);
26117667f90SBarry Smith 
26217667f90SBarry Smith   nzeros = 0;
26317667f90SBarry Smith   ia[0]  = 0;
26417667f90SBarry Smith   for (i=0; i<m; i++) {
26517667f90SBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
26617667f90SBarry Smith     cnt     = 0;
26717667f90SBarry Smith     for (j=0; j<len; j++) {
26817667f90SBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
26917667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27017667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
27117667f90SBarry Smith       }
27217667f90SBarry Smith     }
27317667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27417667f90SBarry Smith     nzeros += cnt;
27517667f90SBarry Smith     ia[i+1] = nzeros;
27617667f90SBarry Smith   }
27717667f90SBarry Smith 
27817667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
27917667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
28058ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28117667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28217667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
28317667f90SBarry Smith 
28417667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
28517667f90SBarry Smith     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
28617667f90SBarry Smith   } else {
28717667f90SBarry Smith     *newmat = B;
28817667f90SBarry Smith   }
28917667f90SBarry Smith   PetscFunctionReturn(0);
29017667f90SBarry Smith }
29117667f90SBarry Smith 
292b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2930b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2943eda8832SBarry Smith        MatGetRow_MPIAdj,
2953eda8832SBarry Smith        MatRestoreRow_MPIAdj,
296b97cf49bSBarry Smith        0,
29797304618SKris Buschelman /* 4*/ 0,
298b97cf49bSBarry Smith        0,
299b97cf49bSBarry Smith        0,
300b97cf49bSBarry Smith        0,
3010b3b1487SBarry Smith        0,
3020b3b1487SBarry Smith        0,
30397304618SKris Buschelman /*10*/ 0,
3040b3b1487SBarry Smith        0,
3050b3b1487SBarry Smith        0,
3060b3b1487SBarry Smith        0,
3070b3b1487SBarry Smith        0,
30897304618SKris Buschelman /*15*/ 0,
3093eda8832SBarry Smith        MatEqual_MPIAdj,
3100b3b1487SBarry Smith        0,
3110b3b1487SBarry Smith        0,
3120b3b1487SBarry Smith        0,
31397304618SKris Buschelman /*20*/ 0,
3140b3b1487SBarry Smith        0,
3153eda8832SBarry Smith        MatSetOption_MPIAdj,
3160b3b1487SBarry Smith        0,
317d519adbfSMatthew Knepley /*24*/ 0,
3180b3b1487SBarry Smith        0,
3190b3b1487SBarry Smith        0,
3200b3b1487SBarry Smith        0,
3210b3b1487SBarry Smith        0,
322d519adbfSMatthew Knepley /*29*/ 0,
3230b3b1487SBarry Smith        0,
324273d9f13SBarry Smith        0,
325273d9f13SBarry Smith        0,
3260b3b1487SBarry Smith        0,
327d519adbfSMatthew Knepley /*34*/ 0,
3280b3b1487SBarry Smith        0,
3290b3b1487SBarry Smith        0,
3300b3b1487SBarry Smith        0,
3310b3b1487SBarry Smith        0,
332d519adbfSMatthew Knepley /*39*/ 0,
3330b3b1487SBarry Smith        0,
3340b3b1487SBarry Smith        0,
3350b3b1487SBarry Smith        0,
3360b3b1487SBarry Smith        0,
337d519adbfSMatthew Knepley /*44*/ 0,
3380b3b1487SBarry Smith        0,
3390b3b1487SBarry Smith        0,
3400b3b1487SBarry Smith        0,
3410b3b1487SBarry Smith        0,
342d519adbfSMatthew Knepley /*49*/ 0,
3436cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
344d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
345b97cf49bSBarry Smith        0,
346b97cf49bSBarry Smith        0,
347d519adbfSMatthew Knepley /*54*/ 0,
348b97cf49bSBarry Smith        0,
3490752156aSBarry Smith        0,
3500752156aSBarry Smith        0,
3510b3b1487SBarry Smith        0,
352d519adbfSMatthew Knepley /*59*/ 0,
353b9b97703SBarry Smith        MatDestroy_MPIAdj,
354b9b97703SBarry Smith        MatView_MPIAdj,
35517667f90SBarry Smith        MatConvertFrom_MPIAdj,
35697304618SKris Buschelman        0,
357d519adbfSMatthew Knepley /*64*/ 0,
35897304618SKris Buschelman        0,
35997304618SKris Buschelman        0,
36097304618SKris Buschelman        0,
36197304618SKris Buschelman        0,
362d519adbfSMatthew Knepley /*69*/ 0,
36397304618SKris Buschelman        0,
36497304618SKris Buschelman        0,
36597304618SKris Buschelman        0,
36697304618SKris Buschelman        0,
367d519adbfSMatthew Knepley /*74*/ 0,
36897304618SKris Buschelman        0,
36997304618SKris Buschelman        0,
37097304618SKris Buschelman        0,
37197304618SKris Buschelman        0,
372d519adbfSMatthew Knepley /*79*/ 0,
37397304618SKris Buschelman        0,
37497304618SKris Buschelman        0,
37597304618SKris Buschelman        0,
376865e5f61SKris Buschelman        0,
377d519adbfSMatthew Knepley /*84*/ 0,
378865e5f61SKris Buschelman        0,
379865e5f61SKris Buschelman        0,
380865e5f61SKris Buschelman        0,
381865e5f61SKris Buschelman        0,
382d519adbfSMatthew Knepley /*89*/ 0,
383865e5f61SKris Buschelman        0,
384865e5f61SKris Buschelman        0,
385865e5f61SKris Buschelman        0,
386865e5f61SKris Buschelman        0,
387d519adbfSMatthew Knepley /*94*/ 0,
388865e5f61SKris Buschelman        0,
389865e5f61SKris Buschelman        0,
390865e5f61SKris Buschelman        0};
391b97cf49bSBarry Smith 
392a23d5eceSKris Buschelman EXTERN_C_BEGIN
393a23d5eceSKris Buschelman #undef __FUNCT__
394a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
3957087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
396a23d5eceSKris Buschelman {
397a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj *)B->data;
398dfbe8321SBarry Smith   PetscErrorCode ierr;
3992515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
400b24ad042SBarry Smith   PetscInt       ii;
401a23d5eceSKris Buschelman #endif
402a23d5eceSKris Buschelman 
403a23d5eceSKris Buschelman   PetscFunctionBegin;
40426283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
40526283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
40626283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
40726283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
40858ec18a5SLisandro Dalcin 
4092515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
410e32f2f54SBarry 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]);
411d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
412a23d5eceSKris Buschelman     if (i[ii] < 0 || i[ii] < i[ii-1]) {
413e32f2f54SBarry Smith       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]);
414a23d5eceSKris Buschelman     }
415a23d5eceSKris Buschelman   }
416d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
417d0f46423SBarry Smith     if (j[ii] < 0 || j[ii] >= B->cmap->N) {
418e32f2f54SBarry Smith       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]);
419a23d5eceSKris Buschelman     }
420a23d5eceSKris Buschelman   }
421a23d5eceSKris Buschelman #endif
42258ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
423a23d5eceSKris Buschelman 
424a23d5eceSKris Buschelman   b->j      = j;
425a23d5eceSKris Buschelman   b->i      = i;
426a23d5eceSKris Buschelman   b->values = values;
427a23d5eceSKris Buschelman 
428d0f46423SBarry Smith   b->nz               = i[B->rmap->n];
429a23d5eceSKris Buschelman   b->diag             = 0;
430a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
431a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
432a23d5eceSKris Buschelman 
433a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
434a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
435a23d5eceSKris Buschelman   PetscFunctionReturn(0);
436a23d5eceSKris Buschelman }
437a23d5eceSKris Buschelman EXTERN_C_END
438b97cf49bSBarry Smith 
4390bad9183SKris Buschelman /*MC
440fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
4410bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
4420bad9183SKris Buschelman 
4430bad9183SKris Buschelman   Level: beginner
4440bad9183SKris Buschelman 
4450bad9183SKris Buschelman .seealso: MatCreateMPIAdj
4460bad9183SKris Buschelman M*/
4470bad9183SKris Buschelman 
448273d9f13SBarry Smith EXTERN_C_BEGIN
4494a2ae208SSatish Balay #undef __FUNCT__
4504a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
4517087cfbeSBarry Smith PetscErrorCode  MatCreate_MPIAdj(Mat B)
452273d9f13SBarry Smith {
453273d9f13SBarry Smith   Mat_MPIAdj     *b;
4546849ba73SBarry Smith   PetscErrorCode ierr;
455273d9f13SBarry Smith 
456273d9f13SBarry Smith   PetscFunctionBegin;
45738f2d2fdSLisandro Dalcin   ierr                = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr);
458b0a32e0cSBarry Smith   B->data             = (void*)b;
459273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
460273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
461273d9f13SBarry Smith 
462a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
463a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
464a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
46517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
466273d9f13SBarry Smith   PetscFunctionReturn(0);
467273d9f13SBarry Smith }
468273d9f13SBarry Smith EXTERN_C_END
469273d9f13SBarry Smith 
4704a2ae208SSatish Balay #undef __FUNCT__
4714a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
472a23d5eceSKris Buschelman /*@C
473a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
474a23d5eceSKris Buschelman 
4753f9fe445SBarry Smith    Logically Collective on MPI_Comm
476a23d5eceSKris Buschelman 
477a23d5eceSKris Buschelman    Input Parameters:
478a23d5eceSKris Buschelman +  A - the matrix
479a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
480a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
481a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
482a23d5eceSKris Buschelman -  values - [optional] edge weights
483a23d5eceSKris Buschelman 
484a23d5eceSKris Buschelman    Level: intermediate
485a23d5eceSKris Buschelman 
486a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
487a23d5eceSKris Buschelman @*/
4887087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
489273d9f13SBarry Smith {
4904ac538c5SBarry Smith   PetscErrorCode ierr;
491273d9f13SBarry Smith 
492273d9f13SBarry Smith   PetscFunctionBegin;
4934ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
494273d9f13SBarry Smith   PetscFunctionReturn(0);
495273d9f13SBarry Smith }
496273d9f13SBarry Smith 
4974a2ae208SSatish Balay #undef __FUNCT__
4984a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
499b97cf49bSBarry Smith /*@C
5003eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
501b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
502b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
503b97cf49bSBarry Smith 
504ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
505ef5ee4d1SLois Curfman McInnes 
506b97cf49bSBarry Smith    Input Parameters:
507c2e958feSBarry Smith +  comm - MPI communicator
5080752156aSBarry Smith .  m - number of local rows
50958ec18a5SLisandro Dalcin .  N - number of global columns
510b97cf49bSBarry Smith .  i - the indices into j for the start of each row
511329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
512ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
513329f5518SBarry Smith -  values -[optional] edge weights
514b97cf49bSBarry Smith 
515b97cf49bSBarry Smith    Output Parameter:
516b97cf49bSBarry Smith .  A - the matrix
517b97cf49bSBarry Smith 
51815091d37SBarry Smith    Level: intermediate
51915091d37SBarry Smith 
5204bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
5214bc6d8c0SBarry Smith    MatSetValues().
522c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
523ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
5241198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
525327e1a73SBarry Smith    Should not include the matrix diagonals.
526b97cf49bSBarry Smith 
527e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
528ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
529ededeb1aSvictorle 
530ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
531b97cf49bSBarry Smith 
532e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
533b97cf49bSBarry Smith @*/
5347087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
535b97cf49bSBarry Smith {
536dfbe8321SBarry Smith   PetscErrorCode ierr;
537b97cf49bSBarry Smith 
538433994e6SBarry Smith   PetscFunctionBegin;
539f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
54058ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
541273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
542273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5433a40ed3dSBarry Smith   PetscFunctionReturn(0);
544b97cf49bSBarry Smith }
545b97cf49bSBarry Smith 
546c2e958feSBarry Smith 
547433994e6SBarry Smith 
548433994e6SBarry Smith 
549433994e6SBarry Smith 
550433994e6SBarry Smith 
551433994e6SBarry Smith 
552