xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision bf0cc55543cd83e035744be2f77202b216d1436e)
1be1d678aSKris Buschelman 
2b97cf49bSBarry Smith /*
3b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
4b97cf49bSBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h>
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) {
23e3c5b3baSBarry Smith     SETERRQ(((PetscObject)A)->comm,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;
492692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5032077d6dSBarry Smith   if (iascii) {
513eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
525cd90555SBarry Smith   } else {
53e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
54b97cf49bSBarry Smith   }
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56b97cf49bSBarry Smith }
57b97cf49bSBarry Smith 
584a2ae208SSatish Balay #undef __FUNCT__
594a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
60dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat)
61b97cf49bSBarry Smith {
623eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
63dfbe8321SBarry Smith   PetscErrorCode ierr;
6494d884c6SBarry Smith 
6594d884c6SBarry Smith   PetscFunctionBegin;
66aa482453SBarry Smith #if defined(PETSC_USE_LOG)
67d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
68b97cf49bSBarry Smith #endif
6905b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
700400557aSBarry Smith   if (a->freeaij) {
7114196391SBarry Smith     if (a->freeaijwithfree) {
7214196391SBarry Smith       if (a->i) free(a->i);
7314196391SBarry Smith       if (a->j) free(a->j);
7414196391SBarry Smith     } else {
75606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
76606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
7705b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
7814196391SBarry Smith     }
790400557aSBarry Smith   }
80*bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
81dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
82901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_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;
12309f38230SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
12409f38230SBarry Smith   ierr = PetscLogObjectMemory(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;
141b24ad042SBarry Smith   PetscInt   *itmp;
142b97cf49bSBarry Smith 
143433994e6SBarry Smith   PetscFunctionBegin;
144d0f46423SBarry Smith   row -= A->rmap->rstart;
1450752156aSBarry Smith 
146e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
147b97cf49bSBarry Smith 
148b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
149b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
150b97cf49bSBarry Smith   if (idx) {
151b97cf49bSBarry Smith     itmp = a->j + a->i[row];
152b97cf49bSBarry Smith     if (*nz) {
153b97cf49bSBarry Smith       *idx = itmp;
154b97cf49bSBarry Smith     }
155b97cf49bSBarry Smith     else *idx = 0;
156b97cf49bSBarry Smith   }
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
158b97cf49bSBarry Smith }
159b97cf49bSBarry Smith 
1604a2ae208SSatish Balay #undef __FUNCT__
1614a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
162b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
163b97cf49bSBarry Smith {
164433994e6SBarry Smith   PetscFunctionBegin;
1653a40ed3dSBarry Smith   PetscFunctionReturn(0);
166b97cf49bSBarry Smith }
167b97cf49bSBarry Smith 
1684a2ae208SSatish Balay #undef __FUNCT__
1694a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
170ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
171b97cf49bSBarry Smith {
1723eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
173dfbe8321SBarry Smith   PetscErrorCode ierr;
174ace3abfcSBarry Smith   PetscBool      flag;
175b97cf49bSBarry Smith 
176433994e6SBarry Smith   PetscFunctionBegin;
177b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
178d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1790f5bd95cSBarry Smith     flag = PETSC_FALSE;
180b97cf49bSBarry Smith   }
181b97cf49bSBarry Smith 
182b97cf49bSBarry Smith   /* if the a->i are the same */
183d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
184b97cf49bSBarry Smith 
185b97cf49bSBarry Smith   /* if a->j are the same */
186b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
187b97cf49bSBarry Smith 
1887adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
1893a40ed3dSBarry Smith   PetscFunctionReturn(0);
190b97cf49bSBarry Smith }
191b97cf49bSBarry Smith 
1924a2ae208SSatish Balay #undef __FUNCT__
1934a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
194ace3abfcSBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool  *done)
1956cd91f64SBarry Smith {
196b24ad042SBarry Smith   PetscInt       i;
1976cd91f64SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
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"
215ace3abfcSBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool  *done)
216d42735a1SBarry Smith {
217b24ad042SBarry Smith   PetscInt   i;
218d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
219d42735a1SBarry Smith 
220d42735a1SBarry Smith   PetscFunctionBegin;
221e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
222e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
223d42735a1SBarry Smith   if (oshift) {
224d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
225d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
226d42735a1SBarry Smith       (*ja)[i]--;
227d42735a1SBarry Smith     }
228d42735a1SBarry Smith   }
2296cd91f64SBarry Smith   PetscFunctionReturn(0);
2306cd91f64SBarry Smith }
231b97cf49bSBarry Smith 
23217667f90SBarry Smith #undef __FUNCT__
23317667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
2347087cfbeSBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
23517667f90SBarry Smith {
23617667f90SBarry Smith   Mat               B;
23717667f90SBarry Smith   PetscErrorCode    ierr;
23817667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
23917667f90SBarry Smith   const PetscInt    *rj;
24017667f90SBarry Smith   const PetscScalar *ra;
24117667f90SBarry Smith   MPI_Comm          comm;
24217667f90SBarry Smith 
24317667f90SBarry Smith   PetscFunctionBegin;
24417667f90SBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
24517667f90SBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
24617667f90SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
24717667f90SBarry Smith 
24817667f90SBarry Smith   /* count the number of nonzeros per row */
24917667f90SBarry Smith   for (i=0; i<m; i++) {
25017667f90SBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
25117667f90SBarry Smith     for (j=0; j<len; j++) {
25217667f90SBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
25317667f90SBarry Smith     }
25417667f90SBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
25517667f90SBarry Smith     nzeros += len;
25617667f90SBarry Smith   }
25717667f90SBarry Smith 
25817667f90SBarry Smith   /* malloc space for nonzeros */
25917667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr);
26017667f90SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
26117667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr);
26217667f90SBarry Smith 
26317667f90SBarry Smith   nzeros = 0;
26417667f90SBarry Smith   ia[0]  = 0;
26517667f90SBarry Smith   for (i=0; i<m; i++) {
26617667f90SBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
26717667f90SBarry Smith     cnt     = 0;
26817667f90SBarry Smith     for (j=0; j<len; j++) {
26917667f90SBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
27017667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27117667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
27217667f90SBarry Smith       }
27317667f90SBarry Smith     }
27417667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27517667f90SBarry Smith     nzeros += cnt;
27617667f90SBarry Smith     ia[i+1] = nzeros;
27717667f90SBarry Smith   }
27817667f90SBarry Smith 
27917667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
28017667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
28158ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28217667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28317667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
28417667f90SBarry Smith 
28517667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
28617667f90SBarry Smith     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
28717667f90SBarry Smith   } else {
28817667f90SBarry Smith     *newmat = B;
28917667f90SBarry Smith   }
29017667f90SBarry Smith   PetscFunctionReturn(0);
29117667f90SBarry Smith }
29217667f90SBarry Smith 
293b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2940b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2953eda8832SBarry Smith        MatGetRow_MPIAdj,
2963eda8832SBarry Smith        MatRestoreRow_MPIAdj,
297b97cf49bSBarry Smith        0,
29897304618SKris Buschelman /* 4*/ 0,
299b97cf49bSBarry Smith        0,
300b97cf49bSBarry Smith        0,
301b97cf49bSBarry Smith        0,
3020b3b1487SBarry Smith        0,
3030b3b1487SBarry Smith        0,
30497304618SKris Buschelman /*10*/ 0,
3050b3b1487SBarry Smith        0,
3060b3b1487SBarry Smith        0,
3070b3b1487SBarry Smith        0,
3080b3b1487SBarry Smith        0,
30997304618SKris Buschelman /*15*/ 0,
3103eda8832SBarry Smith        MatEqual_MPIAdj,
3110b3b1487SBarry Smith        0,
3120b3b1487SBarry Smith        0,
3130b3b1487SBarry Smith        0,
31497304618SKris Buschelman /*20*/ 0,
3150b3b1487SBarry Smith        0,
3163eda8832SBarry Smith        MatSetOption_MPIAdj,
3170b3b1487SBarry Smith        0,
318d519adbfSMatthew Knepley /*24*/ 0,
3190b3b1487SBarry Smith        0,
3200b3b1487SBarry Smith        0,
3210b3b1487SBarry Smith        0,
3220b3b1487SBarry Smith        0,
323d519adbfSMatthew Knepley /*29*/ 0,
3240b3b1487SBarry Smith        0,
325273d9f13SBarry Smith        0,
326273d9f13SBarry Smith        0,
3270b3b1487SBarry Smith        0,
328d519adbfSMatthew Knepley /*34*/ 0,
3290b3b1487SBarry Smith        0,
3300b3b1487SBarry Smith        0,
3310b3b1487SBarry Smith        0,
3320b3b1487SBarry Smith        0,
333d519adbfSMatthew Knepley /*39*/ 0,
3340b3b1487SBarry Smith        0,
3350b3b1487SBarry Smith        0,
3360b3b1487SBarry Smith        0,
3370b3b1487SBarry Smith        0,
338d519adbfSMatthew Knepley /*44*/ 0,
3390b3b1487SBarry Smith        0,
3400b3b1487SBarry Smith        0,
3410b3b1487SBarry Smith        0,
3420b3b1487SBarry Smith        0,
343d519adbfSMatthew Knepley /*49*/ 0,
3446cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
345d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
346b97cf49bSBarry Smith        0,
347b97cf49bSBarry Smith        0,
348d519adbfSMatthew Knepley /*54*/ 0,
349b97cf49bSBarry Smith        0,
3500752156aSBarry Smith        0,
3510752156aSBarry Smith        0,
3520b3b1487SBarry Smith        0,
353d519adbfSMatthew Knepley /*59*/ 0,
354b9b97703SBarry Smith        MatDestroy_MPIAdj,
355b9b97703SBarry Smith        MatView_MPIAdj,
35617667f90SBarry Smith        MatConvertFrom_MPIAdj,
35797304618SKris Buschelman        0,
358d519adbfSMatthew Knepley /*64*/ 0,
35997304618SKris Buschelman        0,
36097304618SKris Buschelman        0,
36197304618SKris Buschelman        0,
36297304618SKris Buschelman        0,
363d519adbfSMatthew Knepley /*69*/ 0,
36497304618SKris Buschelman        0,
36597304618SKris Buschelman        0,
36697304618SKris Buschelman        0,
36797304618SKris Buschelman        0,
368d519adbfSMatthew Knepley /*74*/ 0,
36997304618SKris Buschelman        0,
37097304618SKris Buschelman        0,
37197304618SKris Buschelman        0,
37297304618SKris Buschelman        0,
373d519adbfSMatthew Knepley /*79*/ 0,
37497304618SKris Buschelman        0,
37597304618SKris Buschelman        0,
37697304618SKris Buschelman        0,
377865e5f61SKris Buschelman        0,
378d519adbfSMatthew Knepley /*84*/ 0,
379865e5f61SKris Buschelman        0,
380865e5f61SKris Buschelman        0,
381865e5f61SKris Buschelman        0,
382865e5f61SKris Buschelman        0,
383d519adbfSMatthew Knepley /*89*/ 0,
384865e5f61SKris Buschelman        0,
385865e5f61SKris Buschelman        0,
386865e5f61SKris Buschelman        0,
387865e5f61SKris Buschelman        0,
388d519adbfSMatthew Knepley /*94*/ 0,
389865e5f61SKris Buschelman        0,
390865e5f61SKris Buschelman        0,
391865e5f61SKris Buschelman        0};
392b97cf49bSBarry Smith 
393a23d5eceSKris Buschelman EXTERN_C_BEGIN
394a23d5eceSKris Buschelman #undef __FUNCT__
395a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
3967087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
397a23d5eceSKris Buschelman {
398a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj *)B->data;
399dfbe8321SBarry Smith   PetscErrorCode ierr;
4002515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
401b24ad042SBarry Smith   PetscInt       ii;
402a23d5eceSKris Buschelman #endif
403a23d5eceSKris Buschelman 
404a23d5eceSKris Buschelman   PetscFunctionBegin;
40526283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr);
40626283091SBarry Smith   ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr);
40726283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
40826283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
40958ec18a5SLisandro Dalcin 
4102515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
411e32f2f54SBarry 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]);
412d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
413a23d5eceSKris Buschelman     if (i[ii] < 0 || i[ii] < i[ii-1]) {
414e32f2f54SBarry 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]);
415a23d5eceSKris Buschelman     }
416a23d5eceSKris Buschelman   }
417d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
418d0f46423SBarry Smith     if (j[ii] < 0 || j[ii] >= B->cmap->N) {
419e32f2f54SBarry Smith       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]);
420a23d5eceSKris Buschelman     }
421a23d5eceSKris Buschelman   }
422a23d5eceSKris Buschelman #endif
42358ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
424a23d5eceSKris Buschelman 
425a23d5eceSKris Buschelman   b->j      = j;
426a23d5eceSKris Buschelman   b->i      = i;
427a23d5eceSKris Buschelman   b->values = values;
428a23d5eceSKris Buschelman 
429d0f46423SBarry Smith   b->nz               = i[B->rmap->n];
430a23d5eceSKris Buschelman   b->diag             = 0;
431a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
432a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
433a23d5eceSKris Buschelman 
434a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
435a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
436a23d5eceSKris Buschelman   PetscFunctionReturn(0);
437a23d5eceSKris Buschelman }
438a23d5eceSKris Buschelman EXTERN_C_END
439b97cf49bSBarry Smith 
4400bad9183SKris Buschelman /*MC
441fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
4420bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
4430bad9183SKris Buschelman 
4440bad9183SKris Buschelman   Level: beginner
4450bad9183SKris Buschelman 
4460bad9183SKris Buschelman .seealso: MatCreateMPIAdj
4470bad9183SKris Buschelman M*/
4480bad9183SKris Buschelman 
449273d9f13SBarry Smith EXTERN_C_BEGIN
4504a2ae208SSatish Balay #undef __FUNCT__
4514a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
4527087cfbeSBarry Smith PetscErrorCode  MatCreate_MPIAdj(Mat B)
453273d9f13SBarry Smith {
454273d9f13SBarry Smith   Mat_MPIAdj     *b;
4556849ba73SBarry Smith   PetscErrorCode ierr;
456273d9f13SBarry Smith 
457273d9f13SBarry Smith   PetscFunctionBegin;
45838f2d2fdSLisandro Dalcin   ierr                = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr);
459b0a32e0cSBarry Smith   B->data             = (void*)b;
460273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
461273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
462273d9f13SBarry Smith 
463a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
464a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
465a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
46617667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
467273d9f13SBarry Smith   PetscFunctionReturn(0);
468273d9f13SBarry Smith }
469273d9f13SBarry Smith EXTERN_C_END
470273d9f13SBarry Smith 
4714a2ae208SSatish Balay #undef __FUNCT__
4724a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
473a23d5eceSKris Buschelman /*@C
474a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
475a23d5eceSKris Buschelman 
4763f9fe445SBarry Smith    Logically Collective on MPI_Comm
477a23d5eceSKris Buschelman 
478a23d5eceSKris Buschelman    Input Parameters:
479a23d5eceSKris Buschelman +  A - the matrix
480a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
481a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
482a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
483a23d5eceSKris Buschelman -  values - [optional] edge weights
484a23d5eceSKris Buschelman 
485a23d5eceSKris Buschelman    Level: intermediate
486a23d5eceSKris Buschelman 
487a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
488a23d5eceSKris Buschelman @*/
4897087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
490273d9f13SBarry Smith {
4914ac538c5SBarry Smith   PetscErrorCode ierr;
492273d9f13SBarry Smith 
493273d9f13SBarry Smith   PetscFunctionBegin;
4944ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
495273d9f13SBarry Smith   PetscFunctionReturn(0);
496273d9f13SBarry Smith }
497273d9f13SBarry Smith 
4984a2ae208SSatish Balay #undef __FUNCT__
4994a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
500b97cf49bSBarry Smith /*@C
5013eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
502b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
503b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
504b97cf49bSBarry Smith 
505ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
506ef5ee4d1SLois Curfman McInnes 
507b97cf49bSBarry Smith    Input Parameters:
508c2e958feSBarry Smith +  comm - MPI communicator
5090752156aSBarry Smith .  m - number of local rows
51058ec18a5SLisandro Dalcin .  N - number of global columns
511b97cf49bSBarry Smith .  i - the indices into j for the start of each row
512329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
513ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
514329f5518SBarry Smith -  values -[optional] edge weights
515b97cf49bSBarry Smith 
516b97cf49bSBarry Smith    Output Parameter:
517b97cf49bSBarry Smith .  A - the matrix
518b97cf49bSBarry Smith 
51915091d37SBarry Smith    Level: intermediate
52015091d37SBarry Smith 
5214bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
5224bc6d8c0SBarry Smith    MatSetValues().
523c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
524ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
5251198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
526327e1a73SBarry Smith    Should not include the matrix diagonals.
527b97cf49bSBarry Smith 
528e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
529ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
530ededeb1aSvictorle 
531ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
532b97cf49bSBarry Smith 
533e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
534b97cf49bSBarry Smith @*/
5357087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
536b97cf49bSBarry Smith {
537dfbe8321SBarry Smith   PetscErrorCode ierr;
538b97cf49bSBarry Smith 
539433994e6SBarry Smith   PetscFunctionBegin;
540f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
54158ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
542273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
543273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5443a40ed3dSBarry Smith   PetscFunctionReturn(0);
545b97cf49bSBarry Smith }
546b97cf49bSBarry Smith 
547c2e958feSBarry Smith 
548433994e6SBarry Smith 
549433994e6SBarry Smith 
550433994e6SBarry Smith 
551433994e6SBarry Smith 
552433994e6SBarry Smith 
553