xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 12c028f9705c7f73d22859dd953ec6a759be3484)
1*12c028f9SKris Buschelman /*$Id: mpiadj.c,v 1.61 2001/06/21 21:16:58 bsmith Exp buschelm $*/
2b97cf49bSBarry Smith 
3b97cf49bSBarry Smith /*
4b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
5b97cf49bSBarry Smith */
6e090d566SSatish Balay #include "petscsys.h"
73eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
8b97cf49bSBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
11b0a32e0cSBarry Smith int MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
12b97cf49bSBarry Smith {
133eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
14fb9695e5SSatish Balay   int               ierr,i,j,m = A->m;
15fb9695e5SSatish Balay   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) {
24ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
250752156aSBarry Smith   } else {
26b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
27b97cf49bSBarry Smith     for (i=0; i<m; i++) {
28b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr);
29b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
30b0a32e0cSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr);
310752156aSBarry Smith       }
32b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
33b97cf49bSBarry Smith     }
34b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);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"
42b0a32e0cSBarry Smith int MatView_MPIAdj(Mat A,PetscViewer viewer)
43b97cf49bSBarry Smith {
44b97cf49bSBarry Smith   int        ierr;
456831982aSBarry Smith   PetscTruth isascii;
46b97cf49bSBarry Smith 
47433994e6SBarry Smith   PetscFunctionBegin;
48b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
490f5bd95cSBarry Smith   if (isascii) {
503eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
515cd90555SBarry Smith   } else {
5229bbc08cSBarry Smith     SETERRQ1(1,"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"
593eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat)
60b97cf49bSBarry Smith {
613eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data;
6294d884c6SBarry Smith   int        ierr;
6394d884c6SBarry Smith 
6494d884c6SBarry Smith   PetscFunctionBegin;
65aa482453SBarry Smith #if defined(PETSC_USE_LOG)
66b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
67b97cf49bSBarry Smith #endif
68606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
690400557aSBarry Smith   if (a->freeaij) {
70606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
71606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
721198db28SBarry Smith     if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);}
730400557aSBarry Smith   }
740400557aSBarry Smith   ierr = PetscFree(a->rowners);CHKERRQ(ierr);
751198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
77b97cf49bSBarry Smith }
78b97cf49bSBarry Smith 
794a2ae208SSatish Balay #undef __FUNCT__
804a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
813eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
82b97cf49bSBarry Smith {
833eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
84b97cf49bSBarry Smith 
85433994e6SBarry Smith   PetscFunctionBegin;
86*12c028f9SKris Buschelman   switch (op) {
87*12c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
88b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
89*12c028f9SKris Buschelman     break;
90*12c028f9SKris Buschelman   default:
91b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
92*12c028f9SKris Buschelman     break;
93b97cf49bSBarry Smith   }
943a40ed3dSBarry Smith   PetscFunctionReturn(0);
95b97cf49bSBarry Smith }
96b97cf49bSBarry Smith 
97b97cf49bSBarry Smith 
98b97cf49bSBarry Smith /*
99b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
100b97cf49bSBarry Smith */
101b97cf49bSBarry Smith 
1024a2ae208SSatish Balay #undef __FUNCT__
1034a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
1043eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A)
105b97cf49bSBarry Smith {
1063eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
10782502324SSatish Balay   int        i,j,*diag,m = A->m,ierr;
108b97cf49bSBarry Smith 
109433994e6SBarry Smith   PetscFunctionBegin;
110b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
111b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
112273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
113b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
114b97cf49bSBarry Smith       if (a->j[j] == i) {
115b97cf49bSBarry Smith         diag[i] = j;
116b97cf49bSBarry Smith         break;
117b97cf49bSBarry Smith       }
118b97cf49bSBarry Smith     }
119b97cf49bSBarry Smith   }
120b97cf49bSBarry Smith   a->diag = diag;
1213a40ed3dSBarry Smith   PetscFunctionReturn(0);
122b97cf49bSBarry Smith }
123b97cf49bSBarry Smith 
1244a2ae208SSatish Balay #undef __FUNCT__
1254a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_MPIAdj"
1263eda8832SBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
127b97cf49bSBarry Smith {
1283eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
129433994e6SBarry Smith   PetscFunctionBegin;
130c2e958feSBarry Smith   if (m) *m = a->rstart;
131c2e958feSBarry Smith   if (n) *n = a->rend;
1323a40ed3dSBarry Smith   PetscFunctionReturn(0);
133b97cf49bSBarry Smith }
134b97cf49bSBarry Smith 
1354a2ae208SSatish Balay #undef __FUNCT__
1364a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
1373eda8832SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
138b97cf49bSBarry Smith {
1393eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
140b97cf49bSBarry Smith   int        *itmp;
141b97cf49bSBarry Smith 
142433994e6SBarry Smith   PetscFunctionBegin;
1430752156aSBarry Smith   row -= a->rstart;
1440752156aSBarry Smith 
145273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(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"
1613eda8832SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
162b97cf49bSBarry Smith {
163433994e6SBarry Smith   PetscFunctionBegin;
1643a40ed3dSBarry Smith   PetscFunctionReturn(0);
165b97cf49bSBarry Smith }
166b97cf49bSBarry Smith 
1674a2ae208SSatish Balay #undef __FUNCT__
1684a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj"
1693eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs)
170b97cf49bSBarry Smith {
171433994e6SBarry Smith   PetscFunctionBegin;
172b97cf49bSBarry Smith   *bs = 1;
1733a40ed3dSBarry Smith   PetscFunctionReturn(0);
174b97cf49bSBarry Smith }
175b97cf49bSBarry Smith 
176b97cf49bSBarry Smith 
1774a2ae208SSatish Balay #undef __FUNCT__
1784a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
1793eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
180b97cf49bSBarry Smith {
1813eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
1820f5bd95cSBarry Smith   int         ierr;
1830f5bd95cSBarry Smith   PetscTruth  flag;
184b97cf49bSBarry Smith 
185433994e6SBarry Smith   PetscFunctionBegin;
186273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr);
187273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
188b97cf49bSBarry Smith 
189b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
190273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1910f5bd95cSBarry Smith     flag = PETSC_FALSE;
192b97cf49bSBarry Smith   }
193b97cf49bSBarry Smith 
194b97cf49bSBarry Smith   /* if the a->i are the same */
195273d9f13SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
196b97cf49bSBarry Smith 
197b97cf49bSBarry Smith   /* if a->j are the same */
1980f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
199b97cf49bSBarry Smith 
200ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
2013a40ed3dSBarry Smith   PetscFunctionReturn(0);
202b97cf49bSBarry Smith }
203b97cf49bSBarry Smith 
2044a2ae208SSatish Balay #undef __FUNCT__
2054a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
2066cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
2076cd91f64SBarry Smith {
208d42735a1SBarry Smith   int        ierr,size,i;
2096cd91f64SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
2106cd91f64SBarry Smith 
2116cd91f64SBarry Smith   PetscFunctionBegin;
2126cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2136cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2146cd91f64SBarry Smith   *m    = A->m;
2156cd91f64SBarry Smith   *ia   = a->i;
2166cd91f64SBarry Smith   *ja   = a->j;
2176cd91f64SBarry Smith   *done = PETSC_TRUE;
218d42735a1SBarry Smith   if (oshift) {
219d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
220d42735a1SBarry Smith       (*ja)[i]++;
221d42735a1SBarry Smith     }
222d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
223d42735a1SBarry Smith   }
224d42735a1SBarry Smith   PetscFunctionReturn(0);
225d42735a1SBarry Smith }
226d42735a1SBarry Smith 
2274a2ae208SSatish Balay #undef __FUNCT__
2284a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
229d42735a1SBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
230d42735a1SBarry Smith {
231d42735a1SBarry Smith   int        i;
232d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
233d42735a1SBarry Smith 
234d42735a1SBarry Smith   PetscFunctionBegin;
235c5b3d447SBarry Smith   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
236c5b3d447SBarry Smith   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
237d42735a1SBarry Smith   if (oshift) {
238d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
239d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
240d42735a1SBarry Smith       (*ja)[i]--;
241d42735a1SBarry Smith     }
242d42735a1SBarry Smith   }
2436cd91f64SBarry Smith   PetscFunctionReturn(0);
2446cd91f64SBarry Smith }
245b97cf49bSBarry Smith 
246b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2470b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2483eda8832SBarry Smith        MatGetRow_MPIAdj,
2493eda8832SBarry Smith        MatRestoreRow_MPIAdj,
250b97cf49bSBarry Smith        0,
251b97cf49bSBarry Smith        0,
252b97cf49bSBarry Smith        0,
253b97cf49bSBarry Smith        0,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
2600b3b1487SBarry Smith        0,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2633eda8832SBarry Smith        MatEqual_MPIAdj,
2640b3b1487SBarry Smith        0,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
2670b3b1487SBarry Smith        0,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2703eda8832SBarry Smith        MatSetOption_MPIAdj,
2710b3b1487SBarry Smith        0,
2720b3b1487SBarry Smith        0,
2730b3b1487SBarry Smith        0,
2740b3b1487SBarry Smith        0,
2750b3b1487SBarry Smith        0,
2760b3b1487SBarry Smith        0,
277273d9f13SBarry Smith        0,
278273d9f13SBarry Smith        0,
2793eda8832SBarry Smith        MatGetOwnershipRange_MPIAdj,
2800b3b1487SBarry Smith        0,
2810b3b1487SBarry Smith        0,
2820b3b1487SBarry Smith        0,
2830b3b1487SBarry Smith        0,
2840b3b1487SBarry Smith        0,
2850b3b1487SBarry Smith        0,
2860b3b1487SBarry Smith        0,
2870b3b1487SBarry Smith        0,
2880b3b1487SBarry Smith        0,
2890b3b1487SBarry Smith        0,
2900b3b1487SBarry Smith        0,
2910b3b1487SBarry Smith        0,
2920b3b1487SBarry Smith        0,
2930b3b1487SBarry Smith        0,
2940b3b1487SBarry Smith        0,
2950b3b1487SBarry Smith        0,
2960b3b1487SBarry Smith        0,
2970b3b1487SBarry Smith        0,
298b97cf49bSBarry Smith        0,
2993eda8832SBarry Smith        MatGetBlockSize_MPIAdj,
3006cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
301d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
302b97cf49bSBarry Smith        0,
303b97cf49bSBarry Smith        0,
304b97cf49bSBarry Smith        0,
3050752156aSBarry Smith        0,
3060752156aSBarry Smith        0,
3070b3b1487SBarry Smith        0,
3080b3b1487SBarry Smith        0,
3090b3b1487SBarry Smith        0,
310b9b97703SBarry Smith        MatDestroy_MPIAdj,
311b9b97703SBarry Smith        MatView_MPIAdj,
3120b3b1487SBarry Smith        MatGetMaps_Petsc};
313b97cf49bSBarry Smith 
314b97cf49bSBarry Smith 
315273d9f13SBarry Smith EXTERN_C_BEGIN
3164a2ae208SSatish Balay #undef __FUNCT__
3174a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
318273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B)
319273d9f13SBarry Smith {
320273d9f13SBarry Smith   Mat_MPIAdj *b;
321273d9f13SBarry Smith   int        ii,ierr,size,rank;
322273d9f13SBarry Smith 
323273d9f13SBarry Smith   PetscFunctionBegin;
324273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
325273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
326273d9f13SBarry Smith 
327b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
328b0a32e0cSBarry Smith   B->data             = (void*)b;
329273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
330273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
331273d9f13SBarry Smith   B->factor           = 0;
332273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
333273d9f13SBarry Smith   B->mapping          = 0;
334273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
335273d9f13SBarry Smith 
336273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
337273d9f13SBarry Smith   B->n = B->N;
338273d9f13SBarry Smith 
339273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
340273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
341273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
342273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
343273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
344273d9f13SBarry Smith 
345b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr);
346b0a32e0cSBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
347273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
348273d9f13SBarry Smith   b->rowners[0] = 0;
349273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
350273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
351273d9f13SBarry Smith   }
352273d9f13SBarry Smith   b->rstart = b->rowners[rank];
353273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
354273d9f13SBarry Smith 
355273d9f13SBarry Smith   PetscFunctionReturn(0);
356273d9f13SBarry Smith }
357273d9f13SBarry Smith EXTERN_C_END
358273d9f13SBarry Smith 
3594a2ae208SSatish Balay #undef __FUNCT__
3604a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
361273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
362273d9f13SBarry Smith {
363273d9f13SBarry Smith   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
364273d9f13SBarry Smith   int        ierr;
365273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
366273d9f13SBarry Smith   int        ii;
367273d9f13SBarry Smith #endif
368273d9f13SBarry Smith 
369273d9f13SBarry Smith   PetscFunctionBegin;
370273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
371273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
372273d9f13SBarry Smith   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
373273d9f13SBarry Smith   for (ii=1; ii<B->m; ii++) {
374273d9f13SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) {
375273d9f13SBarry Smith       SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]);
376273d9f13SBarry Smith     }
377273d9f13SBarry Smith   }
378273d9f13SBarry Smith   for (ii=0; ii<i[B->m]; ii++) {
379273d9f13SBarry Smith     if (j[ii] < 0 || j[ii] >= B->N) {
380273d9f13SBarry Smith       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
381273d9f13SBarry Smith     }
382273d9f13SBarry Smith   }
383273d9f13SBarry Smith #endif
384273d9f13SBarry Smith 
385273d9f13SBarry Smith   b->j      = j;
386273d9f13SBarry Smith   b->i      = i;
387273d9f13SBarry Smith   b->values = values;
388273d9f13SBarry Smith 
389273d9f13SBarry Smith   b->nz               = i[B->m];
390273d9f13SBarry Smith   b->diag             = 0;
391273d9f13SBarry Smith   b->symmetric        = PETSC_FALSE;
392273d9f13SBarry Smith   b->freeaij          = PETSC_TRUE;
393273d9f13SBarry Smith 
394273d9f13SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
395273d9f13SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
396273d9f13SBarry Smith   PetscFunctionReturn(0);
397273d9f13SBarry Smith }
398273d9f13SBarry Smith 
3994a2ae208SSatish Balay #undef __FUNCT__
4004a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
401b97cf49bSBarry Smith /*@C
4023eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
403b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
404b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
405b97cf49bSBarry Smith 
406ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
407ef5ee4d1SLois Curfman McInnes 
408b97cf49bSBarry Smith    Input Parameters:
409c2e958feSBarry Smith +  comm - MPI communicator
4100752156aSBarry Smith .  m - number of local rows
411b97cf49bSBarry Smith .  n - number of columns
412b97cf49bSBarry Smith .  i - the indices into j for the start of each row
413329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
414ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
415329f5518SBarry Smith -  values -[optional] edge weights
416b97cf49bSBarry Smith 
417b97cf49bSBarry Smith    Output Parameter:
418b97cf49bSBarry Smith .  A - the matrix
419b97cf49bSBarry Smith 
42015091d37SBarry Smith    Level: intermediate
42115091d37SBarry Smith 
4224bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4234bc6d8c0SBarry Smith    MatSetValues().
424c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
4251198db28SBarry Smith    when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you
4261198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
427327e1a73SBarry Smith    Should not include the matrix diagonals.
428b97cf49bSBarry Smith 
429ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
430b97cf49bSBarry Smith 
43191e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering()
432b97cf49bSBarry Smith @*/
4333eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
434b97cf49bSBarry Smith {
435273d9f13SBarry Smith   int        ierr;
436b97cf49bSBarry Smith 
437433994e6SBarry Smith   PetscFunctionBegin;
438273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
439273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
440273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
4413a40ed3dSBarry Smith   PetscFunctionReturn(0);
442b97cf49bSBarry Smith }
443b97cf49bSBarry Smith 
444273d9f13SBarry Smith EXTERN_C_BEGIN
4454a2ae208SSatish Balay #undef __FUNCT__
4464a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
447273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B)
448c2e958feSBarry Smith {
4494c49b128SBarry Smith   int      i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
450c2e958feSBarry Smith   Scalar   *ra;
451c2e958feSBarry Smith   MPI_Comm comm;
452c2e958feSBarry Smith 
453c2e958feSBarry Smith   PetscFunctionBegin;
454c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
455c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
456c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
457c2e958feSBarry Smith 
458c2e958feSBarry Smith   /* count the number of nonzeros per row */
459c2e958feSBarry Smith   for (i=0; i<m; i++) {
460c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
461c2e958feSBarry Smith     for (j=0; j<len; j++) {
462c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
463c2e958feSBarry Smith     }
464c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
465c2e958feSBarry Smith     nzeros += len;
466c2e958feSBarry Smith   }
467c2e958feSBarry Smith 
468c2e958feSBarry Smith   /* malloc space for nonzeros */
469b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr);
470b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr);
471b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr);
472c2e958feSBarry Smith 
473c2e958feSBarry Smith   nzeros = 0;
474c2e958feSBarry Smith   ia[0]  = 0;
475c2e958feSBarry Smith   for (i=0; i<m; i++) {
476c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
477c2e958feSBarry Smith     cnt     = 0;
478c2e958feSBarry Smith     for (j=0; j<len; j++) {
479c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
4804c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
481c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
482c2e958feSBarry Smith       }
483c2e958feSBarry Smith     }
484c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
485c2e958feSBarry Smith     nzeros += cnt;
486c2e958feSBarry Smith     ia[i+1] = nzeros;
487c2e958feSBarry Smith   }
488c2e958feSBarry Smith 
489c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
4903eda8832SBarry Smith   ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr);
491c2e958feSBarry Smith 
492c2e958feSBarry Smith   PetscFunctionReturn(0);
493c2e958feSBarry Smith }
494273d9f13SBarry Smith EXTERN_C_END
495433994e6SBarry Smith 
496433994e6SBarry Smith 
497433994e6SBarry Smith 
498433994e6SBarry Smith 
499433994e6SBarry Smith 
500