xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 3a7fca6b988f1faea172e46abdf950477dbfaf76)
1*3a7fca6bSBarry Smith /*$Id: mpiadj.c,v 1.60 2001/03/23 23:22:17 balay 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 
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;
19*3a7fca6bSBarry 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;
86b97cf49bSBarry Smith   if (op == MAT_STRUCTURALLY_SYMMETRIC) {
87b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
88b97cf49bSBarry Smith   } else {
89b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
90b97cf49bSBarry Smith   }
913a40ed3dSBarry Smith   PetscFunctionReturn(0);
92b97cf49bSBarry Smith }
93b97cf49bSBarry Smith 
94b97cf49bSBarry Smith 
95b97cf49bSBarry Smith /*
96b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
97b97cf49bSBarry Smith */
98b97cf49bSBarry Smith 
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
1013eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A)
102b97cf49bSBarry Smith {
1033eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
10482502324SSatish Balay   int        i,j,*diag,m = A->m,ierr;
105b97cf49bSBarry Smith 
106433994e6SBarry Smith   PetscFunctionBegin;
107b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
108b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
109273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
110b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
111b97cf49bSBarry Smith       if (a->j[j] == i) {
112b97cf49bSBarry Smith         diag[i] = j;
113b97cf49bSBarry Smith         break;
114b97cf49bSBarry Smith       }
115b97cf49bSBarry Smith     }
116b97cf49bSBarry Smith   }
117b97cf49bSBarry Smith   a->diag = diag;
1183a40ed3dSBarry Smith   PetscFunctionReturn(0);
119b97cf49bSBarry Smith }
120b97cf49bSBarry Smith 
1214a2ae208SSatish Balay #undef __FUNCT__
1224a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_MPIAdj"
1233eda8832SBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
124b97cf49bSBarry Smith {
1253eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
126433994e6SBarry Smith   PetscFunctionBegin;
127c2e958feSBarry Smith   if (m) *m = a->rstart;
128c2e958feSBarry Smith   if (n) *n = a->rend;
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
130b97cf49bSBarry Smith }
131b97cf49bSBarry Smith 
1324a2ae208SSatish Balay #undef __FUNCT__
1334a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
1343eda8832SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
135b97cf49bSBarry Smith {
1363eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
137b97cf49bSBarry Smith   int        *itmp;
138b97cf49bSBarry Smith 
139433994e6SBarry Smith   PetscFunctionBegin;
1400752156aSBarry Smith   row -= a->rstart;
1410752156aSBarry Smith 
142273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
143b97cf49bSBarry Smith 
144b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
145b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
146b97cf49bSBarry Smith   if (idx) {
147b97cf49bSBarry Smith     itmp = a->j + a->i[row];
148b97cf49bSBarry Smith     if (*nz) {
149b97cf49bSBarry Smith       *idx = itmp;
150b97cf49bSBarry Smith     }
151b97cf49bSBarry Smith     else *idx = 0;
152b97cf49bSBarry Smith   }
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154b97cf49bSBarry Smith }
155b97cf49bSBarry Smith 
1564a2ae208SSatish Balay #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
1583eda8832SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
159b97cf49bSBarry Smith {
160433994e6SBarry Smith   PetscFunctionBegin;
1613a40ed3dSBarry Smith   PetscFunctionReturn(0);
162b97cf49bSBarry Smith }
163b97cf49bSBarry Smith 
1644a2ae208SSatish Balay #undef __FUNCT__
1654a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj"
1663eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs)
167b97cf49bSBarry Smith {
168433994e6SBarry Smith   PetscFunctionBegin;
169b97cf49bSBarry Smith   *bs = 1;
1703a40ed3dSBarry Smith   PetscFunctionReturn(0);
171b97cf49bSBarry Smith }
172b97cf49bSBarry Smith 
173b97cf49bSBarry Smith 
1744a2ae208SSatish Balay #undef __FUNCT__
1754a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
1763eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
177b97cf49bSBarry Smith {
1783eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
1790f5bd95cSBarry Smith   int         ierr;
1800f5bd95cSBarry Smith   PetscTruth  flag;
181b97cf49bSBarry Smith 
182433994e6SBarry Smith   PetscFunctionBegin;
183273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr);
184273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
185b97cf49bSBarry Smith 
186b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
187273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1880f5bd95cSBarry Smith     flag = PETSC_FALSE;
189b97cf49bSBarry Smith   }
190b97cf49bSBarry Smith 
191b97cf49bSBarry Smith   /* if the a->i are the same */
192273d9f13SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
193b97cf49bSBarry Smith 
194b97cf49bSBarry Smith   /* if a->j are the same */
1950f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
196b97cf49bSBarry Smith 
197ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
199b97cf49bSBarry Smith }
200b97cf49bSBarry Smith 
2014a2ae208SSatish Balay #undef __FUNCT__
2024a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
2036cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
2046cd91f64SBarry Smith {
205d42735a1SBarry Smith   int        ierr,size,i;
2066cd91f64SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
2076cd91f64SBarry Smith 
2086cd91f64SBarry Smith   PetscFunctionBegin;
2096cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2106cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2116cd91f64SBarry Smith   *m    = A->m;
2126cd91f64SBarry Smith   *ia   = a->i;
2136cd91f64SBarry Smith   *ja   = a->j;
2146cd91f64SBarry Smith   *done = PETSC_TRUE;
215d42735a1SBarry Smith   if (oshift) {
216d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
217d42735a1SBarry Smith       (*ja)[i]++;
218d42735a1SBarry Smith     }
219d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
220d42735a1SBarry Smith   }
221d42735a1SBarry Smith   PetscFunctionReturn(0);
222d42735a1SBarry Smith }
223d42735a1SBarry Smith 
2244a2ae208SSatish Balay #undef __FUNCT__
2254a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
226d42735a1SBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
227d42735a1SBarry Smith {
228d42735a1SBarry Smith   int        i;
229d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
230d42735a1SBarry Smith 
231d42735a1SBarry Smith   PetscFunctionBegin;
232c5b3d447SBarry Smith   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
233c5b3d447SBarry Smith   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
234d42735a1SBarry Smith   if (oshift) {
235d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
236d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
237d42735a1SBarry Smith       (*ja)[i]--;
238d42735a1SBarry Smith     }
239d42735a1SBarry Smith   }
2406cd91f64SBarry Smith   PetscFunctionReturn(0);
2416cd91f64SBarry Smith }
242b97cf49bSBarry Smith 
243b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2440b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2453eda8832SBarry Smith        MatGetRow_MPIAdj,
2463eda8832SBarry Smith        MatRestoreRow_MPIAdj,
247b97cf49bSBarry Smith        0,
248b97cf49bSBarry Smith        0,
249b97cf49bSBarry Smith        0,
250b97cf49bSBarry 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,
2590b3b1487SBarry Smith        0,
2603eda8832SBarry Smith        MatEqual_MPIAdj,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
2640b3b1487SBarry Smith        0,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
2673eda8832SBarry Smith        MatSetOption_MPIAdj,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2700b3b1487SBarry Smith        0,
2710b3b1487SBarry Smith        0,
2720b3b1487SBarry Smith        0,
2730b3b1487SBarry Smith        0,
274273d9f13SBarry Smith        0,
275273d9f13SBarry Smith        0,
2763eda8832SBarry Smith        MatGetOwnershipRange_MPIAdj,
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,
2940b3b1487SBarry Smith        0,
295b97cf49bSBarry Smith        0,
2963eda8832SBarry Smith        MatGetBlockSize_MPIAdj,
2976cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
298d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
299b97cf49bSBarry Smith        0,
300b97cf49bSBarry Smith        0,
301b97cf49bSBarry Smith        0,
3020752156aSBarry Smith        0,
3030752156aSBarry Smith        0,
3040b3b1487SBarry Smith        0,
3050b3b1487SBarry Smith        0,
3060b3b1487SBarry Smith        0,
307b9b97703SBarry Smith        MatDestroy_MPIAdj,
308b9b97703SBarry Smith        MatView_MPIAdj,
3090b3b1487SBarry Smith        MatGetMaps_Petsc};
310b97cf49bSBarry Smith 
311b97cf49bSBarry Smith 
312273d9f13SBarry Smith EXTERN_C_BEGIN
3134a2ae208SSatish Balay #undef __FUNCT__
3144a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
315273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B)
316273d9f13SBarry Smith {
317273d9f13SBarry Smith   Mat_MPIAdj *b;
318273d9f13SBarry Smith   int        ii,ierr,size,rank;
319273d9f13SBarry Smith 
320273d9f13SBarry Smith   PetscFunctionBegin;
321273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
322273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
323273d9f13SBarry Smith 
324b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
325b0a32e0cSBarry Smith   B->data             = (void*)b;
326273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
327273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
328273d9f13SBarry Smith   B->factor           = 0;
329273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
330273d9f13SBarry Smith   B->mapping          = 0;
331273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
332273d9f13SBarry Smith 
333273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
334273d9f13SBarry Smith   B->n = B->N;
335273d9f13SBarry Smith 
336273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
337273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
338273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
339273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
340273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
341273d9f13SBarry Smith 
342b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr);
343b0a32e0cSBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
344273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
345273d9f13SBarry Smith   b->rowners[0] = 0;
346273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
347273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
348273d9f13SBarry Smith   }
349273d9f13SBarry Smith   b->rstart = b->rowners[rank];
350273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
351273d9f13SBarry Smith 
352273d9f13SBarry Smith   PetscFunctionReturn(0);
353273d9f13SBarry Smith }
354273d9f13SBarry Smith EXTERN_C_END
355273d9f13SBarry Smith 
3564a2ae208SSatish Balay #undef __FUNCT__
3574a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
358273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
359273d9f13SBarry Smith {
360273d9f13SBarry Smith   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
361273d9f13SBarry Smith   int        ierr;
362273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
363273d9f13SBarry Smith   int        ii;
364273d9f13SBarry Smith #endif
365273d9f13SBarry Smith 
366273d9f13SBarry Smith   PetscFunctionBegin;
367273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
368273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
369273d9f13SBarry Smith   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
370273d9f13SBarry Smith   for (ii=1; ii<B->m; ii++) {
371273d9f13SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) {
372273d9f13SBarry Smith       SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]);
373273d9f13SBarry Smith     }
374273d9f13SBarry Smith   }
375273d9f13SBarry Smith   for (ii=0; ii<i[B->m]; ii++) {
376273d9f13SBarry Smith     if (j[ii] < 0 || j[ii] >= B->N) {
377273d9f13SBarry Smith       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
378273d9f13SBarry Smith     }
379273d9f13SBarry Smith   }
380273d9f13SBarry Smith #endif
381273d9f13SBarry Smith 
382273d9f13SBarry Smith   b->j      = j;
383273d9f13SBarry Smith   b->i      = i;
384273d9f13SBarry Smith   b->values = values;
385273d9f13SBarry Smith 
386273d9f13SBarry Smith   b->nz               = i[B->m];
387273d9f13SBarry Smith   b->diag             = 0;
388273d9f13SBarry Smith   b->symmetric        = PETSC_FALSE;
389273d9f13SBarry Smith   b->freeaij          = PETSC_TRUE;
390273d9f13SBarry Smith 
391273d9f13SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
392273d9f13SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
393273d9f13SBarry Smith   PetscFunctionReturn(0);
394273d9f13SBarry Smith }
395273d9f13SBarry Smith 
3964a2ae208SSatish Balay #undef __FUNCT__
3974a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
398b97cf49bSBarry Smith /*@C
3993eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
400b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
401b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
402b97cf49bSBarry Smith 
403ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
404ef5ee4d1SLois Curfman McInnes 
405b97cf49bSBarry Smith    Input Parameters:
406c2e958feSBarry Smith +  comm - MPI communicator
4070752156aSBarry Smith .  m - number of local rows
408b97cf49bSBarry Smith .  n - number of columns
409b97cf49bSBarry Smith .  i - the indices into j for the start of each row
410329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
411ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
412329f5518SBarry Smith -  values -[optional] edge weights
413b97cf49bSBarry Smith 
414b97cf49bSBarry Smith    Output Parameter:
415b97cf49bSBarry Smith .  A - the matrix
416b97cf49bSBarry Smith 
41715091d37SBarry Smith    Level: intermediate
41815091d37SBarry Smith 
4194bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4204bc6d8c0SBarry Smith    MatSetValues().
421c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
4221198db28SBarry Smith    when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you
4231198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
424327e1a73SBarry Smith    Should not include the matrix diagonals.
425b97cf49bSBarry Smith 
426ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
427b97cf49bSBarry Smith 
42891e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering()
429b97cf49bSBarry Smith @*/
4303eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
431b97cf49bSBarry Smith {
432273d9f13SBarry Smith   int        ierr;
433b97cf49bSBarry Smith 
434433994e6SBarry Smith   PetscFunctionBegin;
435273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
436273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
437273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
439b97cf49bSBarry Smith }
440b97cf49bSBarry Smith 
441273d9f13SBarry Smith EXTERN_C_BEGIN
4424a2ae208SSatish Balay #undef __FUNCT__
4434a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
444273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B)
445c2e958feSBarry Smith {
4464c49b128SBarry Smith   int      i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
447c2e958feSBarry Smith   Scalar   *ra;
448c2e958feSBarry Smith   MPI_Comm comm;
449c2e958feSBarry Smith 
450c2e958feSBarry Smith   PetscFunctionBegin;
451c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
452c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
453c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
454c2e958feSBarry Smith 
455c2e958feSBarry Smith   /* count the number of nonzeros per row */
456c2e958feSBarry Smith   for (i=0; i<m; i++) {
457c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
458c2e958feSBarry Smith     for (j=0; j<len; j++) {
459c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
460c2e958feSBarry Smith     }
461c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
462c2e958feSBarry Smith     nzeros += len;
463c2e958feSBarry Smith   }
464c2e958feSBarry Smith 
465c2e958feSBarry Smith   /* malloc space for nonzeros */
466b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr);
467b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr);
468b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr);
469c2e958feSBarry Smith 
470c2e958feSBarry Smith   nzeros = 0;
471c2e958feSBarry Smith   ia[0]  = 0;
472c2e958feSBarry Smith   for (i=0; i<m; i++) {
473c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
474c2e958feSBarry Smith     cnt     = 0;
475c2e958feSBarry Smith     for (j=0; j<len; j++) {
476c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
4774c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
478c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
479c2e958feSBarry Smith       }
480c2e958feSBarry Smith     }
481c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
482c2e958feSBarry Smith     nzeros += cnt;
483c2e958feSBarry Smith     ia[i+1] = nzeros;
484c2e958feSBarry Smith   }
485c2e958feSBarry Smith 
486c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
4873eda8832SBarry Smith   ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr);
488c2e958feSBarry Smith 
489c2e958feSBarry Smith   PetscFunctionReturn(0);
490c2e958feSBarry Smith }
491273d9f13SBarry Smith EXTERN_C_END
492433994e6SBarry Smith 
493433994e6SBarry Smith 
494433994e6SBarry Smith 
495433994e6SBarry Smith 
496433994e6SBarry Smith 
497