xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 325e03ae252f1b9b8bd112fd0e09178d7d92063d)
173f4d377SMatthew Knepley /*$Id: mpiadj.c,v 1.66 2001/08/07 03:02:59 balay Exp $*/
2b97cf49bSBarry Smith 
3b97cf49bSBarry Smith /*
4b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
5b97cf49bSBarry Smith */
63eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
7*325e03aeSBarry Smith #include "petscsys.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;
8612c028f9SKris Buschelman   switch (op) {
8712c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
88b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
8912c028f9SKris Buschelman     break;
90d03495bdSKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
9112c028f9SKris Buschelman   default:
92b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
9312c028f9SKris Buschelman     break;
94b97cf49bSBarry Smith   }
953a40ed3dSBarry Smith   PetscFunctionReturn(0);
96b97cf49bSBarry Smith }
97b97cf49bSBarry Smith 
98b97cf49bSBarry Smith 
99b97cf49bSBarry Smith /*
100b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
101b97cf49bSBarry Smith */
102b97cf49bSBarry Smith 
1034a2ae208SSatish Balay #undef __FUNCT__
1044a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
1053eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A)
106b97cf49bSBarry Smith {
1073eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
10882502324SSatish Balay   int        i,j,*diag,m = A->m,ierr;
109b97cf49bSBarry Smith 
110433994e6SBarry Smith   PetscFunctionBegin;
111b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
112b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
113273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
114b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
115b97cf49bSBarry Smith       if (a->j[j] == i) {
116b97cf49bSBarry Smith         diag[i] = j;
117b97cf49bSBarry Smith         break;
118b97cf49bSBarry Smith       }
119b97cf49bSBarry Smith     }
120b97cf49bSBarry Smith   }
121b97cf49bSBarry Smith   a->diag = diag;
1223a40ed3dSBarry Smith   PetscFunctionReturn(0);
123b97cf49bSBarry Smith }
124b97cf49bSBarry Smith 
1254a2ae208SSatish Balay #undef __FUNCT__
1264a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
12787828ca2SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
128b97cf49bSBarry Smith {
1293eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
130b97cf49bSBarry Smith   int        *itmp;
131b97cf49bSBarry Smith 
132433994e6SBarry Smith   PetscFunctionBegin;
1330752156aSBarry Smith   row -= a->rstart;
1340752156aSBarry Smith 
135273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
136b97cf49bSBarry Smith 
137b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
138b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
139b97cf49bSBarry Smith   if (idx) {
140b97cf49bSBarry Smith     itmp = a->j + a->i[row];
141b97cf49bSBarry Smith     if (*nz) {
142b97cf49bSBarry Smith       *idx = itmp;
143b97cf49bSBarry Smith     }
144b97cf49bSBarry Smith     else *idx = 0;
145b97cf49bSBarry Smith   }
1463a40ed3dSBarry Smith   PetscFunctionReturn(0);
147b97cf49bSBarry Smith }
148b97cf49bSBarry Smith 
1494a2ae208SSatish Balay #undef __FUNCT__
1504a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
15187828ca2SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
152b97cf49bSBarry Smith {
153433994e6SBarry Smith   PetscFunctionBegin;
1543a40ed3dSBarry Smith   PetscFunctionReturn(0);
155b97cf49bSBarry Smith }
156b97cf49bSBarry Smith 
1574a2ae208SSatish Balay #undef __FUNCT__
1584a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj"
1593eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs)
160b97cf49bSBarry Smith {
161433994e6SBarry Smith   PetscFunctionBegin;
162b97cf49bSBarry Smith   *bs = 1;
1633a40ed3dSBarry Smith   PetscFunctionReturn(0);
164b97cf49bSBarry Smith }
165b97cf49bSBarry Smith 
166b97cf49bSBarry Smith 
1674a2ae208SSatish Balay #undef __FUNCT__
1684a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
1693eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
170b97cf49bSBarry Smith {
1713eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
1720f5bd95cSBarry Smith   int         ierr;
1730f5bd95cSBarry Smith   PetscTruth  flag;
174b97cf49bSBarry Smith 
175433994e6SBarry Smith   PetscFunctionBegin;
176273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr);
177273d9f13SBarry Smith   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
178b97cf49bSBarry Smith 
179b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
180273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1810f5bd95cSBarry Smith     flag = PETSC_FALSE;
182b97cf49bSBarry Smith   }
183b97cf49bSBarry Smith 
184b97cf49bSBarry Smith   /* if the a->i are the same */
185273d9f13SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
186b97cf49bSBarry Smith 
187b97cf49bSBarry Smith   /* if a->j are the same */
1880f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
189b97cf49bSBarry Smith 
190ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1913a40ed3dSBarry Smith   PetscFunctionReturn(0);
192b97cf49bSBarry Smith }
193b97cf49bSBarry Smith 
1944a2ae208SSatish Balay #undef __FUNCT__
1954a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
1966cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
1976cd91f64SBarry Smith {
198d42735a1SBarry Smith   int        ierr,size,i;
1996cd91f64SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
2006cd91f64SBarry Smith 
2016cd91f64SBarry Smith   PetscFunctionBegin;
2026cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2036cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2046cd91f64SBarry Smith   *m    = A->m;
2056cd91f64SBarry Smith   *ia   = a->i;
2066cd91f64SBarry Smith   *ja   = a->j;
2076cd91f64SBarry Smith   *done = PETSC_TRUE;
208d42735a1SBarry Smith   if (oshift) {
209d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
210d42735a1SBarry Smith       (*ja)[i]++;
211d42735a1SBarry Smith     }
212d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
213d42735a1SBarry Smith   }
214d42735a1SBarry Smith   PetscFunctionReturn(0);
215d42735a1SBarry Smith }
216d42735a1SBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
219d42735a1SBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
220d42735a1SBarry Smith {
221d42735a1SBarry Smith   int        i;
222d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
223d42735a1SBarry Smith 
224d42735a1SBarry Smith   PetscFunctionBegin;
225c5b3d447SBarry Smith   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
226c5b3d447SBarry Smith   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
227d42735a1SBarry Smith   if (oshift) {
228d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
229d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
230d42735a1SBarry Smith       (*ja)[i]--;
231d42735a1SBarry Smith     }
232d42735a1SBarry Smith   }
2336cd91f64SBarry Smith   PetscFunctionReturn(0);
2346cd91f64SBarry Smith }
235b97cf49bSBarry Smith 
236b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2370b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2383eda8832SBarry Smith        MatGetRow_MPIAdj,
2393eda8832SBarry Smith        MatRestoreRow_MPIAdj,
240b97cf49bSBarry Smith        0,
241b97cf49bSBarry Smith        0,
242b97cf49bSBarry Smith        0,
243b97cf49bSBarry Smith        0,
2440b3b1487SBarry Smith        0,
2450b3b1487SBarry Smith        0,
2460b3b1487SBarry Smith        0,
2470b3b1487SBarry Smith        0,
2480b3b1487SBarry Smith        0,
2490b3b1487SBarry Smith        0,
2500b3b1487SBarry Smith        0,
2510b3b1487SBarry Smith        0,
2520b3b1487SBarry Smith        0,
2533eda8832SBarry Smith        MatEqual_MPIAdj,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
2603eda8832SBarry Smith        MatSetOption_MPIAdj,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
2640b3b1487SBarry Smith        0,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
267273d9f13SBarry Smith        0,
268273d9f13SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2700b3b1487SBarry Smith        0,
2710b3b1487SBarry Smith        0,
2720b3b1487SBarry Smith        0,
2730b3b1487SBarry Smith        0,
2740b3b1487SBarry Smith        0,
2750b3b1487SBarry Smith        0,
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,
2873eda8832SBarry Smith        MatGetBlockSize_MPIAdj,
2886cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
289d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
290b97cf49bSBarry Smith        0,
291b97cf49bSBarry Smith        0,
292b97cf49bSBarry Smith        0,
2930752156aSBarry Smith        0,
2940752156aSBarry Smith        0,
2950b3b1487SBarry Smith        0,
2960b3b1487SBarry Smith        0,
2970b3b1487SBarry Smith        0,
298b9b97703SBarry Smith        MatDestroy_MPIAdj,
299b9b97703SBarry Smith        MatView_MPIAdj,
3008a124369SBarry Smith        MatGetPetscMaps_Petsc};
301b97cf49bSBarry Smith 
302b97cf49bSBarry Smith 
303273d9f13SBarry Smith EXTERN_C_BEGIN
3044a2ae208SSatish Balay #undef __FUNCT__
3054a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
306273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B)
307273d9f13SBarry Smith {
308273d9f13SBarry Smith   Mat_MPIAdj *b;
309273d9f13SBarry Smith   int        ii,ierr,size,rank;
310273d9f13SBarry Smith 
311273d9f13SBarry Smith   PetscFunctionBegin;
312273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
313273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
314273d9f13SBarry Smith 
315b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
316b0a32e0cSBarry Smith   B->data             = (void*)b;
317273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
318273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
319273d9f13SBarry Smith   B->factor           = 0;
320273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
321273d9f13SBarry Smith   B->mapping          = 0;
322273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
323273d9f13SBarry Smith 
324273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
325273d9f13SBarry Smith   B->n = B->N;
326273d9f13SBarry Smith 
327273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
328273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
3298a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
330273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
3318a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
332273d9f13SBarry Smith 
333b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr);
334b0a32e0cSBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
335273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
336273d9f13SBarry Smith   b->rowners[0] = 0;
337273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
338273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
339273d9f13SBarry Smith   }
340273d9f13SBarry Smith   b->rstart = b->rowners[rank];
341273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
342273d9f13SBarry Smith 
343273d9f13SBarry Smith   PetscFunctionReturn(0);
344273d9f13SBarry Smith }
345273d9f13SBarry Smith EXTERN_C_END
346273d9f13SBarry Smith 
3474a2ae208SSatish Balay #undef __FUNCT__
3484a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
349273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
350273d9f13SBarry Smith {
351273d9f13SBarry Smith   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
352273d9f13SBarry Smith   int        ierr;
353273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
354273d9f13SBarry Smith   int        ii;
355273d9f13SBarry Smith #endif
356273d9f13SBarry Smith 
357273d9f13SBarry Smith   PetscFunctionBegin;
358273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
359273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
360273d9f13SBarry Smith   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
361273d9f13SBarry Smith   for (ii=1; ii<B->m; ii++) {
362273d9f13SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) {
363273d9f13SBarry Smith       SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]);
364273d9f13SBarry Smith     }
365273d9f13SBarry Smith   }
366273d9f13SBarry Smith   for (ii=0; ii<i[B->m]; ii++) {
367273d9f13SBarry Smith     if (j[ii] < 0 || j[ii] >= B->N) {
368273d9f13SBarry Smith       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
369273d9f13SBarry Smith     }
370273d9f13SBarry Smith   }
371273d9f13SBarry Smith #endif
372273d9f13SBarry Smith 
373273d9f13SBarry Smith   b->j      = j;
374273d9f13SBarry Smith   b->i      = i;
375273d9f13SBarry Smith   b->values = values;
376273d9f13SBarry Smith 
377273d9f13SBarry Smith   b->nz               = i[B->m];
378273d9f13SBarry Smith   b->diag             = 0;
379273d9f13SBarry Smith   b->symmetric        = PETSC_FALSE;
380273d9f13SBarry Smith   b->freeaij          = PETSC_TRUE;
381273d9f13SBarry Smith 
382273d9f13SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
383273d9f13SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
384273d9f13SBarry Smith   PetscFunctionReturn(0);
385273d9f13SBarry Smith }
386273d9f13SBarry Smith 
3874a2ae208SSatish Balay #undef __FUNCT__
3884a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
389b97cf49bSBarry Smith /*@C
3903eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
391b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
392b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
393b97cf49bSBarry Smith 
394ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
395ef5ee4d1SLois Curfman McInnes 
396b97cf49bSBarry Smith    Input Parameters:
397c2e958feSBarry Smith +  comm - MPI communicator
3980752156aSBarry Smith .  m - number of local rows
399b97cf49bSBarry Smith .  n - number of columns
400b97cf49bSBarry Smith .  i - the indices into j for the start of each row
401329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
402ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
403329f5518SBarry Smith -  values -[optional] edge weights
404b97cf49bSBarry Smith 
405b97cf49bSBarry Smith    Output Parameter:
406b97cf49bSBarry Smith .  A - the matrix
407b97cf49bSBarry Smith 
40815091d37SBarry Smith    Level: intermediate
40915091d37SBarry Smith 
4104bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4114bc6d8c0SBarry Smith    MatSetValues().
412c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
4131198db28SBarry Smith    when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you
4141198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
415327e1a73SBarry Smith    Should not include the matrix diagonals.
416b97cf49bSBarry Smith 
417ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
418b97cf49bSBarry Smith 
41991e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering()
420b97cf49bSBarry Smith @*/
4213eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
422b97cf49bSBarry Smith {
423273d9f13SBarry Smith   int        ierr;
424b97cf49bSBarry Smith 
425433994e6SBarry Smith   PetscFunctionBegin;
426273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
427273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
428273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
4293a40ed3dSBarry Smith   PetscFunctionReturn(0);
430b97cf49bSBarry Smith }
431b97cf49bSBarry Smith 
432273d9f13SBarry Smith EXTERN_C_BEGIN
4334a2ae208SSatish Balay #undef __FUNCT__
4344a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
435273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B)
436c2e958feSBarry Smith {
4374c49b128SBarry Smith   int          i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
43887828ca2SBarry Smith   PetscScalar  *ra;
439c2e958feSBarry Smith   MPI_Comm     comm;
440c2e958feSBarry Smith 
441c2e958feSBarry Smith   PetscFunctionBegin;
442c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
443c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
444c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
445c2e958feSBarry Smith 
446c2e958feSBarry Smith   /* count the number of nonzeros per row */
447c2e958feSBarry Smith   for (i=0; i<m; i++) {
448c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
449c2e958feSBarry Smith     for (j=0; j<len; j++) {
450c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
451c2e958feSBarry Smith     }
452c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
453c2e958feSBarry Smith     nzeros += len;
454c2e958feSBarry Smith   }
455c2e958feSBarry Smith 
456c2e958feSBarry Smith   /* malloc space for nonzeros */
457b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr);
458b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr);
459b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr);
460c2e958feSBarry Smith 
461c2e958feSBarry Smith   nzeros = 0;
462c2e958feSBarry Smith   ia[0]  = 0;
463c2e958feSBarry Smith   for (i=0; i<m; i++) {
464c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
465c2e958feSBarry Smith     cnt     = 0;
466c2e958feSBarry Smith     for (j=0; j<len; j++) {
467c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
4684c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
469c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
470c2e958feSBarry Smith       }
471c2e958feSBarry Smith     }
472c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
473c2e958feSBarry Smith     nzeros += cnt;
474c2e958feSBarry Smith     ia[i+1] = nzeros;
475c2e958feSBarry Smith   }
476c2e958feSBarry Smith 
477c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
4783eda8832SBarry Smith   ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr);
479c2e958feSBarry Smith 
480c2e958feSBarry Smith   PetscFunctionReturn(0);
481c2e958feSBarry Smith }
482273d9f13SBarry Smith EXTERN_C_END
483433994e6SBarry Smith 
484433994e6SBarry Smith 
485433994e6SBarry Smith 
486433994e6SBarry Smith 
487433994e6SBarry Smith 
488