xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision b0a32e0c6855ee6a6cd3495fa7da12ea9885bc5d)
1*b0a32e0cSBarry Smith /*$Id: mpiadj.c,v 1.54 2000/11/08 15:15:20 bsmith Exp bsmith $*/
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 
9b97cf49bSBarry Smith #undef __FUNC__
10*b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIAdj_ASCII"
11*b0a32e0cSBarry Smith int MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
12b97cf49bSBarry Smith {
133eda8832SBarry Smith   Mat_MPIAdj  *a = (Mat_MPIAdj*)A->data;
14273d9f13SBarry Smith   int         ierr,i,j,m = A->m, format;
15b97cf49bSBarry Smith   char        *outputname;
16b97cf49bSBarry Smith 
17433994e6SBarry Smith   PetscFunctionBegin;
18*b0a32e0cSBarry Smith   ierr = PetscViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
19*b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
20*b0a32e0cSBarry Smith   if (format == PETSC_VIEWER_FORMAT_ASCII_INFO) {
213a40ed3dSBarry Smith     PetscFunctionReturn(0);
22*b0a32e0cSBarry Smith   } else if (format == PETSC_VIEWER_FORMAT_ASCII_MATLAB) {
23ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
240752156aSBarry Smith   } else {
25*b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
26b97cf49bSBarry Smith     for (i=0; i<m; i++) {
27*b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr);
28b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
29*b0a32e0cSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr);
300752156aSBarry Smith       }
31*b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
32b97cf49bSBarry Smith     }
33*b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
34b97cf49bSBarry Smith   }
35*b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
37b97cf49bSBarry Smith }
38b97cf49bSBarry Smith 
39b97cf49bSBarry Smith #undef __FUNC__
40*b0a32e0cSBarry Smith #define __FUNC__ "MatView_MPIAdj"
41*b0a32e0cSBarry Smith int MatView_MPIAdj(Mat A,PetscViewer viewer)
42b97cf49bSBarry Smith {
43b97cf49bSBarry Smith   int        ierr;
446831982aSBarry Smith   PetscTruth isascii;
45b97cf49bSBarry Smith 
46433994e6SBarry Smith   PetscFunctionBegin;
47*b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
480f5bd95cSBarry Smith   if (isascii) {
493eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
505cd90555SBarry Smith   } else {
5129bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
52b97cf49bSBarry Smith   }
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54b97cf49bSBarry Smith }
55b97cf49bSBarry Smith 
56b97cf49bSBarry Smith #undef __FUNC__
57*b0a32e0cSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj"
583eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat)
59b97cf49bSBarry Smith {
603eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data;
6194d884c6SBarry Smith   int        ierr;
6294d884c6SBarry Smith 
6394d884c6SBarry Smith   PetscFunctionBegin;
64aa482453SBarry Smith #if defined(PETSC_USE_LOG)
65*b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
66b97cf49bSBarry Smith #endif
67606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
680400557aSBarry Smith   if (a->freeaij) {
69606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
70606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
711198db28SBarry Smith     if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);}
720400557aSBarry Smith   }
730400557aSBarry Smith   ierr = PetscFree(a->rowners);CHKERRQ(ierr);
741198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
753a40ed3dSBarry Smith   PetscFunctionReturn(0);
76b97cf49bSBarry Smith }
77b97cf49bSBarry Smith 
78b97cf49bSBarry Smith #undef __FUNC__
79*b0a32e0cSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj"
803eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
81b97cf49bSBarry Smith {
823eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
83b97cf49bSBarry Smith 
84433994e6SBarry Smith   PetscFunctionBegin;
85b97cf49bSBarry Smith   if (op == MAT_STRUCTURALLY_SYMMETRIC) {
86b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
87b97cf49bSBarry Smith   } else {
88*b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
89b97cf49bSBarry Smith   }
903a40ed3dSBarry Smith   PetscFunctionReturn(0);
91b97cf49bSBarry Smith }
92b97cf49bSBarry Smith 
93b97cf49bSBarry Smith 
94b97cf49bSBarry Smith /*
95b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
96b97cf49bSBarry Smith */
97b97cf49bSBarry Smith 
98b97cf49bSBarry Smith #undef __FUNC__
99*b0a32e0cSBarry Smith #define __FUNC__ "MatMarkDiagonal_MPIAdj"
1003eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A)
101b97cf49bSBarry Smith {
1023eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
103273d9f13SBarry Smith   int        i,j,*diag,m = A->m;
104b97cf49bSBarry Smith 
105433994e6SBarry Smith   PetscFunctionBegin;
106*b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&  diag );CHKERRQ(ierr);
107*b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
108273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
109b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
110b97cf49bSBarry Smith       if (a->j[j] == i) {
111b97cf49bSBarry Smith         diag[i] = j;
112b97cf49bSBarry Smith         break;
113b97cf49bSBarry Smith       }
114b97cf49bSBarry Smith     }
115b97cf49bSBarry Smith   }
116b97cf49bSBarry Smith   a->diag = diag;
1173a40ed3dSBarry Smith   PetscFunctionReturn(0);
118b97cf49bSBarry Smith }
119b97cf49bSBarry Smith 
120b97cf49bSBarry Smith #undef __FUNC__
121*b0a32e0cSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj"
1223eda8832SBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
123b97cf49bSBarry Smith {
1243eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
125433994e6SBarry Smith   PetscFunctionBegin;
126c2e958feSBarry Smith   if (m) *m = a->rstart;
127c2e958feSBarry Smith   if (n) *n = a->rend;
1283a40ed3dSBarry Smith   PetscFunctionReturn(0);
129b97cf49bSBarry Smith }
130b97cf49bSBarry Smith 
131b97cf49bSBarry Smith #undef __FUNC__
132*b0a32e0cSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj"
1333eda8832SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
134b97cf49bSBarry Smith {
1353eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
136b97cf49bSBarry Smith   int        *itmp;
137b97cf49bSBarry Smith 
138433994e6SBarry Smith   PetscFunctionBegin;
1390752156aSBarry Smith   row -= a->rstart;
1400752156aSBarry Smith 
141273d9f13SBarry Smith   if (row < 0 || row >= A->m) 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 
155b97cf49bSBarry Smith #undef __FUNC__
156*b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj"
1573eda8832SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
158b97cf49bSBarry Smith {
159433994e6SBarry Smith   PetscFunctionBegin;
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161b97cf49bSBarry Smith }
162b97cf49bSBarry Smith 
163b97cf49bSBarry Smith #undef __FUNC__
164*b0a32e0cSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj"
1653eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs)
166b97cf49bSBarry Smith {
167433994e6SBarry Smith   PetscFunctionBegin;
168b97cf49bSBarry Smith   *bs = 1;
1693a40ed3dSBarry Smith   PetscFunctionReturn(0);
170b97cf49bSBarry Smith }
171b97cf49bSBarry Smith 
172b97cf49bSBarry Smith 
173b97cf49bSBarry Smith #undef __FUNC__
174*b0a32e0cSBarry Smith #define __FUNC__ "MatEqual_MPIAdj"
1753eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
176b97cf49bSBarry Smith {
1773eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
1780f5bd95cSBarry Smith   int         ierr;
1790f5bd95cSBarry Smith   PetscTruth  flag;
180b97cf49bSBarry Smith 
181433994e6SBarry Smith   PetscFunctionBegin;
182273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr);
183273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
184b97cf49bSBarry Smith 
185b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
186273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1870f5bd95cSBarry Smith     flag = PETSC_FALSE;
188b97cf49bSBarry Smith   }
189b97cf49bSBarry Smith 
190b97cf49bSBarry Smith   /* if the a->i are the same */
191273d9f13SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
192b97cf49bSBarry Smith 
193b97cf49bSBarry Smith   /* if a->j are the same */
1940f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
195b97cf49bSBarry Smith 
196ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1973a40ed3dSBarry Smith   PetscFunctionReturn(0);
198b97cf49bSBarry Smith }
199b97cf49bSBarry Smith 
2006cd91f64SBarry Smith #undef __FUNC__
201*b0a32e0cSBarry Smith #define __FUNC__ "MatGetRowIJ_MPIAdj"
2026cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
2036cd91f64SBarry Smith {
204d42735a1SBarry Smith   int        ierr,size,i;
2056cd91f64SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
2066cd91f64SBarry Smith 
2076cd91f64SBarry Smith   PetscFunctionBegin;
2086cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2096cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2106cd91f64SBarry Smith   *m    = A->m;
2116cd91f64SBarry Smith   *ia   = a->i;
2126cd91f64SBarry Smith   *ja   = a->j;
2136cd91f64SBarry Smith   *done = PETSC_TRUE;
214d42735a1SBarry Smith   if (oshift) {
215d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
216d42735a1SBarry Smith       (*ja)[i]++;
217d42735a1SBarry Smith     }
218d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
219d42735a1SBarry Smith   }
220d42735a1SBarry Smith   PetscFunctionReturn(0);
221d42735a1SBarry Smith }
222d42735a1SBarry Smith 
223d42735a1SBarry Smith #undef __FUNC__
224*b0a32e0cSBarry Smith #define __FUNC__ "MatRestoreRowIJ_MPIAdj"
225d42735a1SBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
226d42735a1SBarry Smith {
227d42735a1SBarry Smith   int        i;
228d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
229d42735a1SBarry Smith 
230d42735a1SBarry Smith   PetscFunctionBegin;
231c5b3d447SBarry Smith   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
232c5b3d447SBarry Smith   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
233d42735a1SBarry Smith   if (oshift) {
234d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
235d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
236d42735a1SBarry Smith       (*ja)[i]--;
237d42735a1SBarry Smith     }
238d42735a1SBarry Smith   }
2396cd91f64SBarry Smith   PetscFunctionReturn(0);
2406cd91f64SBarry Smith }
241b97cf49bSBarry Smith 
242b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2430b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2443eda8832SBarry Smith        MatGetRow_MPIAdj,
2453eda8832SBarry Smith        MatRestoreRow_MPIAdj,
246b97cf49bSBarry Smith        0,
247b97cf49bSBarry Smith        0,
248b97cf49bSBarry Smith        0,
249b97cf49bSBarry Smith        0,
2500b3b1487SBarry Smith        0,
2510b3b1487SBarry Smith        0,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        0,
2593eda8832SBarry Smith        MatEqual_MPIAdj,
2600b3b1487SBarry Smith        0,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
2640b3b1487SBarry Smith        0,
2650b3b1487SBarry Smith        0,
2663eda8832SBarry Smith        MatSetOption_MPIAdj,
2670b3b1487SBarry Smith        0,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2700b3b1487SBarry Smith        0,
2710b3b1487SBarry Smith        0,
2720b3b1487SBarry Smith        0,
273273d9f13SBarry Smith        0,
274273d9f13SBarry Smith        0,
2753eda8832SBarry Smith        MatGetOwnershipRange_MPIAdj,
2760b3b1487SBarry Smith        0,
2770b3b1487SBarry Smith        0,
2780b3b1487SBarry Smith        0,
2790b3b1487SBarry Smith        0,
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,
294b97cf49bSBarry Smith        0,
2953eda8832SBarry Smith        MatGetBlockSize_MPIAdj,
2966cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
297d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
298b97cf49bSBarry Smith        0,
299b97cf49bSBarry Smith        0,
300b97cf49bSBarry Smith        0,
3010752156aSBarry Smith        0,
3020752156aSBarry Smith        0,
3030b3b1487SBarry Smith        0,
3040b3b1487SBarry Smith        0,
3050b3b1487SBarry Smith        0,
306b9b97703SBarry Smith        MatDestroy_MPIAdj,
307b9b97703SBarry Smith        MatView_MPIAdj,
3080b3b1487SBarry Smith        MatGetMaps_Petsc};
309b97cf49bSBarry Smith 
310b97cf49bSBarry Smith 
311273d9f13SBarry Smith EXTERN_C_BEGIN
312273d9f13SBarry Smith #undef __FUNC__
313*b0a32e0cSBarry Smith #define __FUNC__ "MatCreate_MPIAdj"
314273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B)
315273d9f13SBarry Smith {
316273d9f13SBarry Smith   Mat_MPIAdj *b;
317273d9f13SBarry Smith   int        ii,ierr,size,rank;
318273d9f13SBarry Smith 
319273d9f13SBarry Smith   PetscFunctionBegin;
320273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
321273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
322273d9f13SBarry Smith 
323*b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
324*b0a32e0cSBarry Smith   B->data             = (void*)b;
325273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
326273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
327273d9f13SBarry Smith   B->factor           = 0;
328273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
329273d9f13SBarry Smith   B->mapping          = 0;
330273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
331273d9f13SBarry Smith 
332273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
333273d9f13SBarry Smith   B->n = B->N;
334273d9f13SBarry Smith 
335273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
336273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
337273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
338273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
339273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
340273d9f13SBarry Smith 
341*b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&  b->rowners );CHKERRQ(ierr);
342*b0a32e0cSBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
343273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
344273d9f13SBarry Smith   b->rowners[0] = 0;
345273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
346273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
347273d9f13SBarry Smith   }
348273d9f13SBarry Smith   b->rstart = b->rowners[rank];
349273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
350273d9f13SBarry Smith 
351273d9f13SBarry Smith   PetscFunctionReturn(0);
352273d9f13SBarry Smith }
353273d9f13SBarry Smith EXTERN_C_END
354273d9f13SBarry Smith 
355273d9f13SBarry Smith #undef __FUNC__
356*b0a32e0cSBarry Smith #define __FUNC__ "MatMPIAdjSetPreallocation"
357273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
358273d9f13SBarry Smith {
359273d9f13SBarry Smith   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
360273d9f13SBarry Smith   int        ierr;
361273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
362273d9f13SBarry Smith   int        ii;
363273d9f13SBarry Smith #endif
364273d9f13SBarry Smith 
365273d9f13SBarry Smith   PetscFunctionBegin;
366273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
367273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
368273d9f13SBarry Smith   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
369273d9f13SBarry Smith   for (ii=1; ii<B->m; ii++) {
370273d9f13SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) {
371273d9f13SBarry Smith       SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]);
372273d9f13SBarry Smith     }
373273d9f13SBarry Smith   }
374273d9f13SBarry Smith   for (ii=0; ii<i[B->m]; ii++) {
375273d9f13SBarry Smith     if (j[ii] < 0 || j[ii] >= B->N) {
376273d9f13SBarry Smith       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
377273d9f13SBarry Smith     }
378273d9f13SBarry Smith   }
379273d9f13SBarry Smith #endif
380273d9f13SBarry Smith 
381273d9f13SBarry Smith   b->j      = j;
382273d9f13SBarry Smith   b->i      = i;
383273d9f13SBarry Smith   b->values = values;
384273d9f13SBarry Smith 
385273d9f13SBarry Smith   b->nz               = i[B->m];
386273d9f13SBarry Smith   b->diag             = 0;
387273d9f13SBarry Smith   b->symmetric        = PETSC_FALSE;
388273d9f13SBarry Smith   b->freeaij          = PETSC_TRUE;
389273d9f13SBarry Smith 
390273d9f13SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
391273d9f13SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
392273d9f13SBarry Smith   PetscFunctionReturn(0);
393273d9f13SBarry Smith }
394273d9f13SBarry Smith 
395b97cf49bSBarry Smith #undef __FUNC__
396*b0a32e0cSBarry Smith #define __FUNC__ "MatCreateMPIAdj"
397b97cf49bSBarry Smith /*@C
3983eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
399b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
400b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
401b97cf49bSBarry Smith 
402ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
403ef5ee4d1SLois Curfman McInnes 
404b97cf49bSBarry Smith    Input Parameters:
405c2e958feSBarry Smith +  comm - MPI communicator
4060752156aSBarry Smith .  m - number of local rows
407b97cf49bSBarry Smith .  n - number of columns
408b97cf49bSBarry Smith .  i - the indices into j for the start of each row
409329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
410ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
411329f5518SBarry Smith -  values -[optional] edge weights
412b97cf49bSBarry Smith 
413b97cf49bSBarry Smith    Output Parameter:
414b97cf49bSBarry Smith .  A - the matrix
415b97cf49bSBarry Smith 
41615091d37SBarry Smith    Level: intermediate
41715091d37SBarry Smith 
4184bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4194bc6d8c0SBarry Smith    MatSetValues().
420c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
4211198db28SBarry Smith    when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you
4221198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
423327e1a73SBarry Smith    Should not include the matrix diagonals.
424b97cf49bSBarry Smith 
425ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
426b97cf49bSBarry Smith 
42791e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering()
428b97cf49bSBarry Smith @*/
4293eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
430b97cf49bSBarry Smith {
431273d9f13SBarry Smith   int        ierr;
432b97cf49bSBarry Smith 
433433994e6SBarry Smith   PetscFunctionBegin;
434273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
435273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
436273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
4373a40ed3dSBarry Smith   PetscFunctionReturn(0);
438b97cf49bSBarry Smith }
439b97cf49bSBarry Smith 
440273d9f13SBarry Smith EXTERN_C_BEGIN
441c2e958feSBarry Smith #undef __FUNC__
442*b0a32e0cSBarry Smith #define __FUNC__ "MatConvertTo_MPIAdj"
443273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B)
444c2e958feSBarry Smith {
4454c49b128SBarry Smith   int      i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
446c2e958feSBarry Smith   Scalar   *ra;
447c2e958feSBarry Smith   MPI_Comm comm;
448c2e958feSBarry Smith 
449c2e958feSBarry Smith   PetscFunctionBegin;
450c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
451c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
452c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
453c2e958feSBarry Smith 
454c2e958feSBarry Smith   /* count the number of nonzeros per row */
455c2e958feSBarry Smith   for (i=0; i<m; i++) {
456c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
457c2e958feSBarry Smith     for (j=0; j<len; j++) {
458c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
459c2e958feSBarry Smith     }
460c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
461c2e958feSBarry Smith     nzeros += len;
462c2e958feSBarry Smith   }
463c2e958feSBarry Smith 
464c2e958feSBarry Smith   /* malloc space for nonzeros */
465*b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&  a  );CHKERRQ(ierr);
466*b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&  ia );CHKERRQ(ierr);
467*b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&  ja );CHKERRQ(ierr);
468c2e958feSBarry Smith 
469c2e958feSBarry Smith   nzeros = 0;
470c2e958feSBarry Smith   ia[0]  = 0;
471c2e958feSBarry Smith   for (i=0; i<m; i++) {
472c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
473c2e958feSBarry Smith     cnt     = 0;
474c2e958feSBarry Smith     for (j=0; j<len; j++) {
475c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
4764c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
477c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
478c2e958feSBarry Smith       }
479c2e958feSBarry Smith     }
480c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
481c2e958feSBarry Smith     nzeros += cnt;
482c2e958feSBarry Smith     ia[i+1] = nzeros;
483c2e958feSBarry Smith   }
484c2e958feSBarry Smith 
485c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
4863eda8832SBarry Smith   ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr);
487c2e958feSBarry Smith 
488c2e958feSBarry Smith   PetscFunctionReturn(0);
489c2e958feSBarry Smith }
490273d9f13SBarry Smith EXTERN_C_END
491433994e6SBarry Smith 
492433994e6SBarry Smith 
493433994e6SBarry Smith 
494433994e6SBarry Smith 
495433994e6SBarry Smith 
496