xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision d519adbfaddd20300d0b94bbdb1ee93b9df0dacd)
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"
7325e03aeSBarry Smith #include "petscsys.h"
8b97cf49bSBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
11dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
12b97cf49bSBarry Smith {
133eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
14dfbe8321SBarry Smith   PetscErrorCode    ierr;
15d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n;
162dcb1b2aSMatthew Knepley   const char        *name;
17ce1411ecSBarry Smith   PetscViewerFormat format;
18b97cf49bSBarry Smith 
19433994e6SBarry Smith   PetscFunctionBegin;
203a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
21b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
22fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
233a40ed3dSBarry Smith     PetscFunctionReturn(0);
24fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
25ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
260752156aSBarry Smith   } else {
27b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
28b97cf49bSBarry Smith     for (i=0; i<m; i++) {
29d0f46423SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr);
30b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
3177431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
320752156aSBarry Smith       }
33b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
34b97cf49bSBarry Smith     }
35b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
36b97cf49bSBarry Smith   }
37b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
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;
4632077d6dSBarry Smith   PetscTruth     iascii;
47b97cf49bSBarry Smith 
48433994e6SBarry Smith   PetscFunctionBegin;
4932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
5032077d6dSBarry Smith   if (iascii) {
513eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
525cd90555SBarry Smith   } else {
5379a5c55eSBarry Smith     SETERRQ1(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) {
71606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
72606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
7305b42c5fSBarry Smith     ierr = PetscFree(a->values);CHKERRQ(ierr);
740400557aSBarry Smith   }
751198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
76dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
77901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
783a40ed3dSBarry Smith   PetscFunctionReturn(0);
79b97cf49bSBarry Smith }
80b97cf49bSBarry Smith 
814a2ae208SSatish Balay #undef __FUNCT__
824a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
834e0d8c25SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscTruth flg)
84b97cf49bSBarry Smith {
853eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
8663ba0a88SBarry Smith   PetscErrorCode ierr;
87b97cf49bSBarry Smith 
88433994e6SBarry Smith   PetscFunctionBegin;
8912c028f9SKris Buschelman   switch (op) {
9077e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9112c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
929a4540c5SBarry Smith   case MAT_HERMITIAN:
934e0d8c25SBarry Smith     a->symmetric = flg;
949a4540c5SBarry Smith     break;
959a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
969a4540c5SBarry Smith     break;
9712c028f9SKris Buschelman   default:
98290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
9912c028f9SKris Buschelman     break;
100b97cf49bSBarry Smith   }
1013a40ed3dSBarry Smith   PetscFunctionReturn(0);
102b97cf49bSBarry Smith }
103b97cf49bSBarry Smith 
104b97cf49bSBarry Smith 
105b97cf49bSBarry Smith /*
106b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
107b97cf49bSBarry Smith */
108b97cf49bSBarry Smith 
1094a2ae208SSatish Balay #undef __FUNCT__
1104a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
111dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
112b97cf49bSBarry Smith {
1133eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1146849ba73SBarry Smith   PetscErrorCode ierr;
115d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n;
116b97cf49bSBarry Smith 
117433994e6SBarry Smith   PetscFunctionBegin;
11809f38230SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
11909f38230SBarry Smith   ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr);
120d0f46423SBarry Smith   for (i=0; i<A->rmap->n; i++) {
121b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
122b97cf49bSBarry Smith       if (a->j[j] == i) {
12309f38230SBarry Smith         a->diag[i] = j;
124b97cf49bSBarry Smith         break;
125b97cf49bSBarry Smith       }
126b97cf49bSBarry Smith     }
127b97cf49bSBarry Smith   }
1283a40ed3dSBarry Smith   PetscFunctionReturn(0);
129b97cf49bSBarry Smith }
130b97cf49bSBarry Smith 
1314a2ae208SSatish Balay #undef __FUNCT__
1324a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
133b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
134b97cf49bSBarry Smith {
1353eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
136b24ad042SBarry Smith   PetscInt   *itmp;
137b97cf49bSBarry Smith 
138433994e6SBarry Smith   PetscFunctionBegin;
139d0f46423SBarry Smith   row -= A->rmap->rstart;
1400752156aSBarry Smith 
141d0f46423SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
142b97cf49bSBarry Smith 
143b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
144b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
145b97cf49bSBarry Smith   if (idx) {
146b97cf49bSBarry Smith     itmp = a->j + a->i[row];
147b97cf49bSBarry Smith     if (*nz) {
148b97cf49bSBarry Smith       *idx = itmp;
149b97cf49bSBarry Smith     }
150b97cf49bSBarry Smith     else *idx = 0;
151b97cf49bSBarry Smith   }
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
153b97cf49bSBarry Smith }
154b97cf49bSBarry Smith 
1554a2ae208SSatish Balay #undef __FUNCT__
1564a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
157b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
158b97cf49bSBarry Smith {
159433994e6SBarry Smith   PetscFunctionBegin;
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161b97cf49bSBarry Smith }
162b97cf49bSBarry Smith 
1634a2ae208SSatish Balay #undef __FUNCT__
1644a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
165dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
166b97cf49bSBarry Smith {
1673eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
168dfbe8321SBarry Smith   PetscErrorCode ierr;
1690f5bd95cSBarry Smith   PetscTruth     flag;
170b97cf49bSBarry Smith 
171433994e6SBarry Smith   PetscFunctionBegin;
172b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
173d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1740f5bd95cSBarry Smith     flag = PETSC_FALSE;
175b97cf49bSBarry Smith   }
176b97cf49bSBarry Smith 
177b97cf49bSBarry Smith   /* if the a->i are the same */
178d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
179b97cf49bSBarry Smith 
180b97cf49bSBarry Smith   /* if a->j are the same */
181b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
182b97cf49bSBarry Smith 
1837adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
1843a40ed3dSBarry Smith   PetscFunctionReturn(0);
185b97cf49bSBarry Smith }
186b97cf49bSBarry Smith 
1874a2ae208SSatish Balay #undef __FUNCT__
1884a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
1898f7157efSSatish Balay PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
1906cd91f64SBarry Smith {
191dfbe8321SBarry Smith   PetscErrorCode ierr;
192b24ad042SBarry Smith   PetscMPIInt    size;
193b24ad042SBarry Smith   PetscInt       i;
1946cd91f64SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
1956cd91f64SBarry Smith 
1966cd91f64SBarry Smith   PetscFunctionBegin;
1977adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1986cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
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"
2148f7157efSSatish Balay PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
215d42735a1SBarry Smith {
216b24ad042SBarry Smith   PetscInt   i;
217d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
218d42735a1SBarry Smith 
219d42735a1SBarry Smith   PetscFunctionBegin;
220e005ede5SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
221e005ede5SBarry Smith   if (ja && a->j != *ja) SETERRQ(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"
233a313700dSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT 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,
317*d519adbfSMatthew Knepley /*24*/ 0,
3180b3b1487SBarry Smith        0,
3190b3b1487SBarry Smith        0,
3200b3b1487SBarry Smith        0,
3210b3b1487SBarry Smith        0,
322*d519adbfSMatthew Knepley /*29*/ 0,
3230b3b1487SBarry Smith        0,
324273d9f13SBarry Smith        0,
325273d9f13SBarry Smith        0,
3260b3b1487SBarry Smith        0,
327*d519adbfSMatthew Knepley /*34*/ 0,
3280b3b1487SBarry Smith        0,
3290b3b1487SBarry Smith        0,
3300b3b1487SBarry Smith        0,
3310b3b1487SBarry Smith        0,
332*d519adbfSMatthew Knepley /*39*/ 0,
3330b3b1487SBarry Smith        0,
3340b3b1487SBarry Smith        0,
3350b3b1487SBarry Smith        0,
3360b3b1487SBarry Smith        0,
337*d519adbfSMatthew Knepley /*44*/ 0,
3380b3b1487SBarry Smith        0,
3390b3b1487SBarry Smith        0,
3400b3b1487SBarry Smith        0,
3410b3b1487SBarry Smith        0,
342*d519adbfSMatthew Knepley /*49*/ 0,
3436cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
344d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
345b97cf49bSBarry Smith        0,
346b97cf49bSBarry Smith        0,
347*d519adbfSMatthew Knepley /*54*/ 0,
348b97cf49bSBarry Smith        0,
3490752156aSBarry Smith        0,
3500752156aSBarry Smith        0,
3510b3b1487SBarry Smith        0,
352*d519adbfSMatthew Knepley /*59*/ 0,
353b9b97703SBarry Smith        MatDestroy_MPIAdj,
354b9b97703SBarry Smith        MatView_MPIAdj,
35517667f90SBarry Smith        MatConvertFrom_MPIAdj,
35697304618SKris Buschelman        0,
357*d519adbfSMatthew Knepley /*64*/ 0,
35897304618SKris Buschelman        0,
35997304618SKris Buschelman        0,
36097304618SKris Buschelman        0,
36197304618SKris Buschelman        0,
362*d519adbfSMatthew Knepley /*69*/ 0,
36397304618SKris Buschelman        0,
36497304618SKris Buschelman        0,
36597304618SKris Buschelman        0,
36697304618SKris Buschelman        0,
367*d519adbfSMatthew Knepley /*74*/ 0,
36897304618SKris Buschelman        0,
36997304618SKris Buschelman        0,
37097304618SKris Buschelman        0,
37197304618SKris Buschelman        0,
372*d519adbfSMatthew Knepley /*79*/ 0,
37397304618SKris Buschelman        0,
37497304618SKris Buschelman        0,
37597304618SKris Buschelman        0,
376865e5f61SKris Buschelman        0,
377*d519adbfSMatthew Knepley /*84*/ 0,
378865e5f61SKris Buschelman        0,
379865e5f61SKris Buschelman        0,
380865e5f61SKris Buschelman        0,
381865e5f61SKris Buschelman        0,
382*d519adbfSMatthew Knepley /*89*/ 0,
383865e5f61SKris Buschelman        0,
384865e5f61SKris Buschelman        0,
385865e5f61SKris Buschelman        0,
386865e5f61SKris Buschelman        0,
387*d519adbfSMatthew 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"
395be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT 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;
40458ec18a5SLisandro Dalcin   ierr = PetscMapSetBlockSize(B->rmap,1);CHKERRQ(ierr);
40558ec18a5SLisandro Dalcin   ierr = PetscMapSetBlockSize(B->cmap,1);CHKERRQ(ierr);
40658ec18a5SLisandro Dalcin   ierr = PetscMapSetUp(B->rmap);CHKERRQ(ierr);
40758ec18a5SLisandro Dalcin   ierr = PetscMapSetUp(B->cmap);CHKERRQ(ierr);
40858ec18a5SLisandro Dalcin 
4092515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
41077431f27SBarry Smith   if (i[0] != 0) SETERRQ1(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]) {
41377431f27SBarry Smith       SETERRQ4(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) {
41877431f27SBarry Smith       SETERRQ2(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"
451be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIAdj(Mat B)
452273d9f13SBarry Smith {
453273d9f13SBarry Smith   Mat_MPIAdj     *b;
4546849ba73SBarry Smith   PetscErrorCode ierr;
455b24ad042SBarry Smith   PetscMPIInt    size,rank;
456273d9f13SBarry Smith 
457273d9f13SBarry Smith   PetscFunctionBegin;
4587adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr);
4597adad957SLisandro Dalcin   ierr = MPI_Comm_rank(((PetscObject)B)->comm,&rank);CHKERRQ(ierr);
460273d9f13SBarry Smith 
46138f2d2fdSLisandro Dalcin   ierr                = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr);
462b0a32e0cSBarry Smith   B->data             = (void*)b;
463273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
464273d9f13SBarry Smith   B->mapping          = 0;
465273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
466273d9f13SBarry Smith 
467a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
468a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
469a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
47017667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
471273d9f13SBarry Smith   PetscFunctionReturn(0);
472273d9f13SBarry Smith }
473273d9f13SBarry Smith EXTERN_C_END
474273d9f13SBarry Smith 
4754a2ae208SSatish Balay #undef __FUNCT__
4764a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
477a23d5eceSKris Buschelman /*@C
478a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
479a23d5eceSKris Buschelman 
480a23d5eceSKris Buschelman    Collective on MPI_Comm
481a23d5eceSKris Buschelman 
482a23d5eceSKris Buschelman    Input Parameters:
483a23d5eceSKris Buschelman +  A - the matrix
484a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
485a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
486a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
487a23d5eceSKris Buschelman -  values - [optional] edge weights
488a23d5eceSKris Buschelman 
489a23d5eceSKris Buschelman    Level: intermediate
490a23d5eceSKris Buschelman 
491a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
492a23d5eceSKris Buschelman @*/
493be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
494273d9f13SBarry Smith {
495b24ad042SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,PetscInt*);
496273d9f13SBarry Smith 
497273d9f13SBarry Smith   PetscFunctionBegin;
498a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
499a23d5eceSKris Buschelman   if (f) {
500a23d5eceSKris Buschelman     ierr = (*f)(B,i,j,values);CHKERRQ(ierr);
501273d9f13SBarry Smith   }
502273d9f13SBarry Smith   PetscFunctionReturn(0);
503273d9f13SBarry Smith }
504273d9f13SBarry Smith 
5054a2ae208SSatish Balay #undef __FUNCT__
5064a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
507b97cf49bSBarry Smith /*@C
5083eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
509b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
510b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
511b97cf49bSBarry Smith 
512ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
513ef5ee4d1SLois Curfman McInnes 
514b97cf49bSBarry Smith    Input Parameters:
515c2e958feSBarry Smith +  comm - MPI communicator
5160752156aSBarry Smith .  m - number of local rows
51758ec18a5SLisandro Dalcin .  N - number of global columns
518b97cf49bSBarry Smith .  i - the indices into j for the start of each row
519329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
520ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
521329f5518SBarry Smith -  values -[optional] edge weights
522b97cf49bSBarry Smith 
523b97cf49bSBarry Smith    Output Parameter:
524b97cf49bSBarry Smith .  A - the matrix
525b97cf49bSBarry Smith 
52615091d37SBarry Smith    Level: intermediate
52715091d37SBarry Smith 
5284bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
5294bc6d8c0SBarry Smith    MatSetValues().
530c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
531ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
5321198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
533327e1a73SBarry Smith    Should not include the matrix diagonals.
534b97cf49bSBarry Smith 
535e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
536ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
537ededeb1aSvictorle 
538ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
539b97cf49bSBarry Smith 
540e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
541b97cf49bSBarry Smith @*/
54258ec18a5SLisandro Dalcin PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
543b97cf49bSBarry Smith {
544dfbe8321SBarry Smith   PetscErrorCode ierr;
545b97cf49bSBarry Smith 
546433994e6SBarry Smith   PetscFunctionBegin;
547f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
54858ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
549273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
550273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5513a40ed3dSBarry Smith   PetscFunctionReturn(0);
552b97cf49bSBarry Smith }
553b97cf49bSBarry Smith 
554c2e958feSBarry Smith 
555433994e6SBarry Smith 
556433994e6SBarry Smith 
557433994e6SBarry Smith 
558433994e6SBarry Smith 
559433994e6SBarry Smith 
560