xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 1a83f524638a4da8317a6bd80eb6d9a2936d8384)
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) {
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;
49251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((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   }
80bf0cc555SLisandro 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);
839a3665c6SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C","",PETSC_NULL);CHKERRQ(ierr);
843a40ed3dSBarry Smith   PetscFunctionReturn(0);
85b97cf49bSBarry Smith }
86b97cf49bSBarry Smith 
874a2ae208SSatish Balay #undef __FUNCT__
884a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
89ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool  flg)
90b97cf49bSBarry Smith {
913eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
9263ba0a88SBarry Smith   PetscErrorCode ierr;
93b97cf49bSBarry Smith 
94433994e6SBarry Smith   PetscFunctionBegin;
9512c028f9SKris Buschelman   switch (op) {
9677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9712c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
989a4540c5SBarry Smith   case MAT_HERMITIAN:
994e0d8c25SBarry Smith     a->symmetric = flg;
1009a4540c5SBarry Smith     break;
1019a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1029a4540c5SBarry Smith     break;
10312c028f9SKris Buschelman   default:
104290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
10512c028f9SKris Buschelman     break;
106b97cf49bSBarry Smith   }
1073a40ed3dSBarry Smith   PetscFunctionReturn(0);
108b97cf49bSBarry Smith }
109b97cf49bSBarry Smith 
110b97cf49bSBarry Smith 
111b97cf49bSBarry Smith /*
112b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
113b97cf49bSBarry Smith */
114b97cf49bSBarry Smith 
1154a2ae208SSatish Balay #undef __FUNCT__
1164a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
117dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
118b97cf49bSBarry Smith {
1193eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1206849ba73SBarry Smith   PetscErrorCode ierr;
121d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n;
122b97cf49bSBarry Smith 
123433994e6SBarry Smith   PetscFunctionBegin;
12409f38230SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
12509f38230SBarry Smith   ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr);
126d0f46423SBarry Smith   for (i=0; i<A->rmap->n; i++) {
127b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
128b97cf49bSBarry Smith       if (a->j[j] == i) {
12909f38230SBarry Smith         a->diag[i] = j;
130b97cf49bSBarry Smith         break;
131b97cf49bSBarry Smith       }
132b97cf49bSBarry Smith     }
133b97cf49bSBarry Smith   }
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
135b97cf49bSBarry Smith }
136b97cf49bSBarry Smith 
1374a2ae208SSatish Balay #undef __FUNCT__
1384a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
139b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
140b97cf49bSBarry Smith {
1413eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
142b24ad042SBarry Smith   PetscInt   *itmp;
143b97cf49bSBarry Smith 
144433994e6SBarry Smith   PetscFunctionBegin;
145d0f46423SBarry Smith   row -= A->rmap->rstart;
1460752156aSBarry Smith 
147e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
148b97cf49bSBarry Smith 
149b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
150b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
151b97cf49bSBarry Smith   if (idx) {
152b97cf49bSBarry Smith     itmp = a->j + a->i[row];
153b97cf49bSBarry Smith     if (*nz) {
154b97cf49bSBarry Smith       *idx = itmp;
155b97cf49bSBarry Smith     }
156b97cf49bSBarry Smith     else *idx = 0;
157b97cf49bSBarry Smith   }
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
159b97cf49bSBarry Smith }
160b97cf49bSBarry Smith 
1614a2ae208SSatish Balay #undef __FUNCT__
1624a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
163b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
164b97cf49bSBarry Smith {
165433994e6SBarry Smith   PetscFunctionBegin;
1663a40ed3dSBarry Smith   PetscFunctionReturn(0);
167b97cf49bSBarry Smith }
168b97cf49bSBarry Smith 
1694a2ae208SSatish Balay #undef __FUNCT__
1704a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
171ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
172b97cf49bSBarry Smith {
1733eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
174dfbe8321SBarry Smith   PetscErrorCode ierr;
175ace3abfcSBarry Smith   PetscBool      flag;
176b97cf49bSBarry Smith 
177433994e6SBarry Smith   PetscFunctionBegin;
178b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
179d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1800f5bd95cSBarry Smith     flag = PETSC_FALSE;
181b97cf49bSBarry Smith   }
182b97cf49bSBarry Smith 
183b97cf49bSBarry Smith   /* if the a->i are the same */
184d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
185b97cf49bSBarry Smith 
186b97cf49bSBarry Smith   /* if a->j are the same */
187b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
188b97cf49bSBarry Smith 
1897adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
1903a40ed3dSBarry Smith   PetscFunctionReturn(0);
191b97cf49bSBarry Smith }
192b97cf49bSBarry Smith 
1934a2ae208SSatish Balay #undef __FUNCT__
1944a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
195*1a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
1966cd91f64SBarry Smith {
197b24ad042SBarry Smith   PetscInt       i;
1986cd91f64SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
199*1a83f524SJed Brown   PetscInt       **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
2006cd91f64SBarry Smith 
2016cd91f64SBarry Smith   PetscFunctionBegin;
202d0f46423SBarry Smith   *m    = A->rmap->n;
2036cd91f64SBarry Smith   *ia   = a->i;
2046cd91f64SBarry Smith   *ja   = a->j;
2056cd91f64SBarry Smith   *done = PETSC_TRUE;
206d42735a1SBarry Smith   if (oshift) {
207d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
208d42735a1SBarry Smith       (*ja)[i]++;
209d42735a1SBarry Smith     }
210d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
211d42735a1SBarry Smith   }
212d42735a1SBarry Smith   PetscFunctionReturn(0);
213d42735a1SBarry Smith }
214d42735a1SBarry Smith 
2154a2ae208SSatish Balay #undef __FUNCT__
2164a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
217*1a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
218d42735a1SBarry Smith {
219b24ad042SBarry Smith   PetscInt   i;
220d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
221*1a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
222d42735a1SBarry Smith 
223d42735a1SBarry Smith   PetscFunctionBegin;
224e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
225e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
226d42735a1SBarry Smith   if (oshift) {
227d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
228d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
229d42735a1SBarry Smith       (*ja)[i]--;
230d42735a1SBarry Smith     }
231d42735a1SBarry Smith   }
2326cd91f64SBarry Smith   PetscFunctionReturn(0);
2336cd91f64SBarry Smith }
234b97cf49bSBarry Smith 
23517667f90SBarry Smith #undef __FUNCT__
23617667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
23719fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
23817667f90SBarry Smith {
23917667f90SBarry Smith   Mat               B;
24017667f90SBarry Smith   PetscErrorCode    ierr;
24117667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
24217667f90SBarry Smith   const PetscInt    *rj;
24317667f90SBarry Smith   const PetscScalar *ra;
24417667f90SBarry Smith   MPI_Comm          comm;
24517667f90SBarry Smith 
24617667f90SBarry Smith   PetscFunctionBegin;
24717667f90SBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
24817667f90SBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
24917667f90SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
25017667f90SBarry Smith 
25117667f90SBarry Smith   /* count the number of nonzeros per row */
25217667f90SBarry Smith   for (i=0; i<m; i++) {
253f8b03baaSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
254f8b03baaSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,PETSC_NULL,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       a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27017667f90SBarry Smith       ja[nzeros+cnt++] = rj[j];
27117667f90SBarry Smith     }
27217667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27317667f90SBarry Smith     nzeros += cnt;
27417667f90SBarry Smith     ia[i+1] = nzeros;
27517667f90SBarry Smith   }
27617667f90SBarry Smith 
27717667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
27817667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
27958ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28017667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28117667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
28217667f90SBarry Smith 
28317667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
28417667f90SBarry Smith     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
28517667f90SBarry Smith   } else {
28617667f90SBarry Smith     *newmat = B;
28717667f90SBarry Smith   }
28817667f90SBarry Smith   PetscFunctionReturn(0);
28917667f90SBarry Smith }
29017667f90SBarry Smith 
291b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2920b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2933eda8832SBarry Smith        MatGetRow_MPIAdj,
2943eda8832SBarry Smith        MatRestoreRow_MPIAdj,
295b97cf49bSBarry Smith        0,
29697304618SKris Buschelman /* 4*/ 0,
297b97cf49bSBarry Smith        0,
298b97cf49bSBarry Smith        0,
299b97cf49bSBarry Smith        0,
3000b3b1487SBarry Smith        0,
3010b3b1487SBarry Smith        0,
30297304618SKris Buschelman /*10*/ 0,
3030b3b1487SBarry Smith        0,
3040b3b1487SBarry Smith        0,
3050b3b1487SBarry Smith        0,
3060b3b1487SBarry Smith        0,
30797304618SKris Buschelman /*15*/ 0,
3083eda8832SBarry Smith        MatEqual_MPIAdj,
3090b3b1487SBarry Smith        0,
3100b3b1487SBarry Smith        0,
3110b3b1487SBarry Smith        0,
31297304618SKris Buschelman /*20*/ 0,
3130b3b1487SBarry Smith        0,
3143eda8832SBarry Smith        MatSetOption_MPIAdj,
3150b3b1487SBarry Smith        0,
316d519adbfSMatthew Knepley /*24*/ 0,
3170b3b1487SBarry Smith        0,
3180b3b1487SBarry Smith        0,
3190b3b1487SBarry Smith        0,
3200b3b1487SBarry Smith        0,
321d519adbfSMatthew Knepley /*29*/ 0,
3220b3b1487SBarry Smith        0,
323273d9f13SBarry Smith        0,
324273d9f13SBarry Smith        0,
3250b3b1487SBarry Smith        0,
326d519adbfSMatthew Knepley /*34*/ 0,
3270b3b1487SBarry Smith        0,
3280b3b1487SBarry Smith        0,
3290b3b1487SBarry Smith        0,
3300b3b1487SBarry Smith        0,
331d519adbfSMatthew Knepley /*39*/ 0,
3320b3b1487SBarry Smith        0,
3330b3b1487SBarry Smith        0,
3340b3b1487SBarry Smith        0,
3350b3b1487SBarry Smith        0,
336d519adbfSMatthew Knepley /*44*/ 0,
3370b3b1487SBarry Smith        0,
3380b3b1487SBarry Smith        0,
3390b3b1487SBarry Smith        0,
3400b3b1487SBarry Smith        0,
341d519adbfSMatthew Knepley /*49*/ 0,
3426cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
343d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
344b97cf49bSBarry Smith        0,
345b97cf49bSBarry Smith        0,
346d519adbfSMatthew Knepley /*54*/ 0,
347b97cf49bSBarry Smith        0,
3480752156aSBarry Smith        0,
3490752156aSBarry Smith        0,
3500b3b1487SBarry Smith        0,
351d519adbfSMatthew Knepley /*59*/ 0,
352b9b97703SBarry Smith        MatDestroy_MPIAdj,
353b9b97703SBarry Smith        MatView_MPIAdj,
35417667f90SBarry Smith        MatConvertFrom_MPIAdj,
35597304618SKris Buschelman        0,
356d519adbfSMatthew Knepley /*64*/ 0,
35797304618SKris Buschelman        0,
35897304618SKris Buschelman        0,
35997304618SKris Buschelman        0,
36097304618SKris Buschelman        0,
361d519adbfSMatthew Knepley /*69*/ 0,
36297304618SKris Buschelman        0,
36397304618SKris Buschelman        0,
36497304618SKris Buschelman        0,
36597304618SKris Buschelman        0,
366d519adbfSMatthew Knepley /*74*/ 0,
36797304618SKris Buschelman        0,
36897304618SKris Buschelman        0,
36997304618SKris Buschelman        0,
37097304618SKris Buschelman        0,
371d519adbfSMatthew Knepley /*79*/ 0,
37297304618SKris Buschelman        0,
37397304618SKris Buschelman        0,
37497304618SKris Buschelman        0,
375865e5f61SKris Buschelman        0,
376d519adbfSMatthew Knepley /*84*/ 0,
377865e5f61SKris Buschelman        0,
378865e5f61SKris Buschelman        0,
379865e5f61SKris Buschelman        0,
380865e5f61SKris Buschelman        0,
381d519adbfSMatthew Knepley /*89*/ 0,
382865e5f61SKris Buschelman        0,
383865e5f61SKris Buschelman        0,
384865e5f61SKris Buschelman        0,
385865e5f61SKris Buschelman        0,
386d519adbfSMatthew Knepley /*94*/ 0,
387865e5f61SKris Buschelman        0,
388865e5f61SKris Buschelman        0,
389865e5f61SKris Buschelman        0};
390b97cf49bSBarry Smith 
391a23d5eceSKris Buschelman EXTERN_C_BEGIN
392a23d5eceSKris Buschelman #undef __FUNCT__
393a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
3947087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
395a23d5eceSKris Buschelman {
396a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj *)B->data;
397dfbe8321SBarry Smith   PetscErrorCode ierr;
3982515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
399b24ad042SBarry Smith   PetscInt       ii;
400a23d5eceSKris Buschelman #endif
401a23d5eceSKris Buschelman 
402a23d5eceSKris Buschelman   PetscFunctionBegin;
40326283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
40426283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
40558ec18a5SLisandro Dalcin 
4062515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
407e32f2f54SBarry 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]);
408d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
409e02043d6SBarry 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]);
410a23d5eceSKris Buschelman   }
411d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
412e02043d6SBarry 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]);
413a23d5eceSKris Buschelman   }
414a23d5eceSKris Buschelman #endif
41558ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
416a23d5eceSKris Buschelman 
417a23d5eceSKris Buschelman   b->j      = j;
418a23d5eceSKris Buschelman   b->i      = i;
419a23d5eceSKris Buschelman   b->values = values;
420a23d5eceSKris Buschelman 
421d0f46423SBarry Smith   b->nz               = i[B->rmap->n];
422a23d5eceSKris Buschelman   b->diag             = 0;
423a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
424a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
425a23d5eceSKris Buschelman 
426a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
427a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
428a23d5eceSKris Buschelman   PetscFunctionReturn(0);
429a23d5eceSKris Buschelman }
430a23d5eceSKris Buschelman EXTERN_C_END
431b97cf49bSBarry Smith 
4329a3665c6SJed Brown #undef __FUNCT__
4339a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj"
4349a3665c6SJed Brown PETSC_EXTERN_C PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
4359a3665c6SJed Brown {
4369a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
4379a3665c6SJed Brown   PetscErrorCode ierr;
4389a3665c6SJed Brown   const PetscInt *ranges;
4399a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
4409a3665c6SJed Brown   MPI_Group      agroup,bgroup;
4419a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
4429a3665c6SJed Brown 
4439a3665c6SJed Brown   PetscFunctionBegin;
4449a3665c6SJed Brown   *B = PETSC_NULL;
4459a3665c6SJed Brown   acomm = ((PetscObject)A)->comm;
4469a3665c6SJed Brown   ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
4479a3665c6SJed Brown   ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
4489a3665c6SJed Brown   ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
4499a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
4509a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
4519a3665c6SJed Brown   }
4529a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
4539a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
4549a3665c6SJed Brown     *B = A;
4559a3665c6SJed Brown     PetscFunctionReturn(0);
4569a3665c6SJed Brown   }
4579a3665c6SJed Brown 
4589a3665c6SJed Brown   ierr = PetscMalloc(nranks*sizeof(PetscMPIInt),&ranks);CHKERRQ(ierr);
4599a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
4609a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
4619a3665c6SJed Brown   }
4629a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
4639a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
4649a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
4659a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
4669a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
4679a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
4689a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
4699a3665c6SJed Brown     PetscInt   m,N;
4709a3665c6SJed Brown     Mat_MPIAdj *b;
4719a3665c6SJed Brown     ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
4729a3665c6SJed Brown     ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
4739a3665c6SJed Brown     ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
4749a3665c6SJed Brown     b = (Mat_MPIAdj*)(*B)->data;
4759a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
4769a3665c6SJed Brown     ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
4779a3665c6SJed Brown   }
4789a3665c6SJed Brown   PetscFunctionReturn(0);
4799a3665c6SJed Brown }
4809a3665c6SJed Brown 
4819a3665c6SJed Brown #undef __FUNCT__
4829a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat"
4839a3665c6SJed Brown /*@
4849a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
4859a3665c6SJed Brown 
4869a3665c6SJed Brown    Collective
4879a3665c6SJed Brown 
4889a3665c6SJed Brown    Input Arguments:
4899a3665c6SJed Brown .  A - original MPIAdj matrix
4909a3665c6SJed Brown 
4919a3665c6SJed Brown    Output Arguments:
4929a3665c6SJed Brown .  B - matrix on subcommunicator, PETSC_NULL on ranks that owned zero rows of A
4939a3665c6SJed Brown 
4949a3665c6SJed Brown    Level: developer
4959a3665c6SJed Brown 
4969a3665c6SJed Brown    Note:
4979a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
4989a3665c6SJed Brown 
4999a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
5009a3665c6SJed Brown 
5019a3665c6SJed Brown .seealso: MatCreateMPIAdj()
5029a3665c6SJed Brown @*/
5039a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
5049a3665c6SJed Brown {
5059a3665c6SJed Brown   PetscErrorCode ierr;
5069a3665c6SJed Brown 
5079a3665c6SJed Brown   PetscFunctionBegin;
5089a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
5099a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
5109a3665c6SJed Brown   PetscFunctionReturn(0);
5119a3665c6SJed Brown }
5129a3665c6SJed Brown 
5130bad9183SKris Buschelman /*MC
514fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
5150bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
5160bad9183SKris Buschelman 
5170bad9183SKris Buschelman   Level: beginner
5180bad9183SKris Buschelman 
5190bad9183SKris Buschelman .seealso: MatCreateMPIAdj
5200bad9183SKris Buschelman M*/
5210bad9183SKris Buschelman 
522273d9f13SBarry Smith EXTERN_C_BEGIN
5234a2ae208SSatish Balay #undef __FUNCT__
5244a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
5257087cfbeSBarry Smith PetscErrorCode  MatCreate_MPIAdj(Mat B)
526273d9f13SBarry Smith {
527273d9f13SBarry Smith   Mat_MPIAdj     *b;
5286849ba73SBarry Smith   PetscErrorCode ierr;
529273d9f13SBarry Smith 
530273d9f13SBarry Smith   PetscFunctionBegin;
53138f2d2fdSLisandro Dalcin   ierr                = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr);
532b0a32e0cSBarry Smith   B->data             = (void*)b;
533273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
534273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
535273d9f13SBarry Smith 
536a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
537a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
538a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
5399a3665c6SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",
5409a3665c6SJed Brown                                     "MatMPIAdjCreateNonemptySubcommMat_MPIAdj",
5419a3665c6SJed Brown                                      MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
54217667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
543273d9f13SBarry Smith   PetscFunctionReturn(0);
544273d9f13SBarry Smith }
545273d9f13SBarry Smith EXTERN_C_END
546273d9f13SBarry Smith 
5474a2ae208SSatish Balay #undef __FUNCT__
5484a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
549a23d5eceSKris Buschelman /*@C
550a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
551a23d5eceSKris Buschelman 
5523f9fe445SBarry Smith    Logically Collective on MPI_Comm
553a23d5eceSKris Buschelman 
554a23d5eceSKris Buschelman    Input Parameters:
555a23d5eceSKris Buschelman +  A - the matrix
556a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
557a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
558a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
559a23d5eceSKris Buschelman -  values - [optional] edge weights
560a23d5eceSKris Buschelman 
561a23d5eceSKris Buschelman    Level: intermediate
562a23d5eceSKris Buschelman 
563a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
564a23d5eceSKris Buschelman @*/
5657087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
566273d9f13SBarry Smith {
5674ac538c5SBarry Smith   PetscErrorCode ierr;
568273d9f13SBarry Smith 
569273d9f13SBarry Smith   PetscFunctionBegin;
5704ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
571273d9f13SBarry Smith   PetscFunctionReturn(0);
572273d9f13SBarry Smith }
573273d9f13SBarry Smith 
5744a2ae208SSatish Balay #undef __FUNCT__
5754a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
576b97cf49bSBarry Smith /*@C
5773eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
578b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
579b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
580b97cf49bSBarry Smith 
581ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
582ef5ee4d1SLois Curfman McInnes 
583b97cf49bSBarry Smith    Input Parameters:
584c2e958feSBarry Smith +  comm - MPI communicator
5850752156aSBarry Smith .  m - number of local rows
58658ec18a5SLisandro Dalcin .  N - number of global columns
587b97cf49bSBarry Smith .  i - the indices into j for the start of each row
588329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
589ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
590329f5518SBarry Smith -  values -[optional] edge weights
591b97cf49bSBarry Smith 
592b97cf49bSBarry Smith    Output Parameter:
593b97cf49bSBarry Smith .  A - the matrix
594b97cf49bSBarry Smith 
59515091d37SBarry Smith    Level: intermediate
59615091d37SBarry Smith 
5974bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
5984bc6d8c0SBarry Smith    MatSetValues().
599c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
600ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
6011198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
602327e1a73SBarry Smith    Should not include the matrix diagonals.
603b97cf49bSBarry Smith 
604e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
605ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
606ededeb1aSvictorle 
607ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
608b97cf49bSBarry Smith 
609e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
610b97cf49bSBarry Smith @*/
6117087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
612b97cf49bSBarry Smith {
613dfbe8321SBarry Smith   PetscErrorCode ierr;
614b97cf49bSBarry Smith 
615433994e6SBarry Smith   PetscFunctionBegin;
616f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
61758ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
618273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
619273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
621b97cf49bSBarry Smith }
622b97cf49bSBarry Smith 
623c2e958feSBarry Smith 
624433994e6SBarry Smith 
625433994e6SBarry Smith 
626433994e6SBarry Smith 
627433994e6SBarry Smith 
628433994e6SBarry Smith 
629