xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision d42735a103fa0eebdc40ab0086e89c7b7740dbb2)
1*d42735a1SBarry Smith /*$Id: mpiadj.c,v 1.51 2000/11/06 16:27:58 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__
10b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIAdj_ASCII"
11ca44d042SBarry Smith int MatView_MPIAdj_ASCII(Mat A,Viewer 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;
1877ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
19d132466eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
20b97cf49bSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
213a40ed3dSBarry Smith     PetscFunctionReturn(0);
22ffac6cdbSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
23ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
240752156aSBarry Smith   } else {
256831982aSBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
26b97cf49bSBarry Smith     for (i=0; i<m; i++) {
276831982aSBarry Smith       ierr = ViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr);
28b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
296831982aSBarry Smith         ierr = ViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr);
300752156aSBarry Smith       }
316831982aSBarry Smith       ierr = ViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
32b97cf49bSBarry Smith     }
336831982aSBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
34b97cf49bSBarry Smith   }
356831982aSBarry Smith   ierr = ViewerFlush(viewer);CHKERRQ(ierr);
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
37b97cf49bSBarry Smith }
38b97cf49bSBarry Smith 
39b97cf49bSBarry Smith #undef __FUNC__
40b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIAdj"
413eda8832SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer)
42b97cf49bSBarry Smith {
43b97cf49bSBarry Smith   int        ierr;
446831982aSBarry Smith   PetscTruth isascii;
45b97cf49bSBarry Smith 
46433994e6SBarry Smith   PetscFunctionBegin;
476831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&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__
57b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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)
6594d884c6SBarry Smith   PLogObjectState((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__
79b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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 {
883eda8832SBarry Smith     PLogInfo(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__
99b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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;
106b97cf49bSBarry Smith   diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag);
107b97cf49bSBarry Smith   PLogObjectMemory(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__
121b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
132b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
156b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
164b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
174b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"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__
2016cd91f64SBarry Smith #define __FUNC__ /*<a name="MatGetRowIJ_MPIAdj"></a>*/"MatGetRowIJ_MPIAdj"
2026cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
2036cd91f64SBarry Smith {
204*d42735a1SBarry 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;
214*d42735a1SBarry Smith   if (oshift) {
215*d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
216*d42735a1SBarry Smith       (*ja)[i]++;
217*d42735a1SBarry Smith     }
218*d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
219*d42735a1SBarry Smith   }
220*d42735a1SBarry Smith   PetscFunctionReturn(0);
221*d42735a1SBarry Smith }
222*d42735a1SBarry Smith 
223*d42735a1SBarry Smith #undef __FUNC__
224*d42735a1SBarry Smith #define __FUNC__ /*<a name="MatRestoreRowIJ_MPIAdj"></a>*/"MatRestoreRowIJ_MPIAdj"
225*d42735a1SBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
226*d42735a1SBarry Smith {
227*d42735a1SBarry Smith   int        i;
228*d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
229*d42735a1SBarry Smith 
230*d42735a1SBarry Smith   PetscFunctionBegin;
231*d42735a1SBarry Smith   if (oshift) {
232*d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
233*d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
234*d42735a1SBarry Smith       (*ja)[i]--;
235*d42735a1SBarry Smith     }
236*d42735a1SBarry Smith   }
2376cd91f64SBarry Smith   PetscFunctionReturn(0);
2386cd91f64SBarry Smith }
239b97cf49bSBarry Smith 
240b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2410b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2423eda8832SBarry Smith        MatGetRow_MPIAdj,
2433eda8832SBarry Smith        MatRestoreRow_MPIAdj,
244b97cf49bSBarry Smith        0,
245b97cf49bSBarry Smith        0,
246b97cf49bSBarry Smith        0,
247b97cf49bSBarry Smith        0,
2480b3b1487SBarry Smith        0,
2490b3b1487SBarry Smith        0,
2500b3b1487SBarry Smith        0,
2510b3b1487SBarry Smith        0,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
2573eda8832SBarry Smith        MatEqual_MPIAdj,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
2600b3b1487SBarry Smith        0,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
2643eda8832SBarry Smith        MatSetOption_MPIAdj,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
2670b3b1487SBarry Smith        0,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2700b3b1487SBarry Smith        0,
271273d9f13SBarry Smith        0,
272273d9f13SBarry Smith        0,
2733eda8832SBarry Smith        MatGetOwnershipRange_MPIAdj,
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,
2870b3b1487SBarry Smith        0,
2880b3b1487SBarry Smith        0,
2890b3b1487SBarry Smith        0,
2900b3b1487SBarry Smith        0,
2910b3b1487SBarry Smith        0,
292b97cf49bSBarry Smith        0,
2933eda8832SBarry Smith        MatGetBlockSize_MPIAdj,
2946cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
295*d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
296b97cf49bSBarry Smith        0,
297b97cf49bSBarry Smith        0,
298b97cf49bSBarry Smith        0,
2990752156aSBarry Smith        0,
3000752156aSBarry Smith        0,
3010b3b1487SBarry Smith        0,
3020b3b1487SBarry Smith        0,
3030b3b1487SBarry Smith        0,
304b9b97703SBarry Smith        MatDestroy_MPIAdj,
305b9b97703SBarry Smith        MatView_MPIAdj,
3060b3b1487SBarry Smith        MatGetMaps_Petsc};
307b97cf49bSBarry Smith 
308b97cf49bSBarry Smith 
309273d9f13SBarry Smith EXTERN_C_BEGIN
310273d9f13SBarry Smith #undef __FUNC__
311273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_MPIAdj"
312273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B)
313273d9f13SBarry Smith {
314273d9f13SBarry Smith   Mat_MPIAdj *b;
315273d9f13SBarry Smith   int        ii,ierr,size,rank;
316273d9f13SBarry Smith 
317273d9f13SBarry Smith   PetscFunctionBegin;
318273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
319273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
320273d9f13SBarry Smith 
321273d9f13SBarry Smith   B->data             = (void*)(b = PetscNew(Mat_MPIAdj));CHKPTRQ(b);
322273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
323273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
324273d9f13SBarry Smith   B->factor           = 0;
325273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
326273d9f13SBarry Smith   B->mapping          = 0;
327273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
328273d9f13SBarry Smith 
329273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
330273d9f13SBarry Smith   B->n = B->N;
331273d9f13SBarry Smith 
332273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
333273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
334273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
335273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
336273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
337273d9f13SBarry Smith 
338273d9f13SBarry Smith   b->rowners = (int*)PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners);
339273d9f13SBarry Smith   PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
340273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
341273d9f13SBarry Smith   b->rowners[0] = 0;
342273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
343273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
344273d9f13SBarry Smith   }
345273d9f13SBarry Smith   b->rstart = b->rowners[rank];
346273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
347273d9f13SBarry Smith 
348273d9f13SBarry Smith   PetscFunctionReturn(0);
349273d9f13SBarry Smith }
350273d9f13SBarry Smith EXTERN_C_END
351273d9f13SBarry Smith 
352273d9f13SBarry Smith #undef __FUNC__
353273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMPIAdjSetPreallocation"
354273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
355273d9f13SBarry Smith {
356273d9f13SBarry Smith   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
357273d9f13SBarry Smith   int        ierr;
358273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
359273d9f13SBarry Smith   int        ii;
360273d9f13SBarry Smith #endif
361273d9f13SBarry Smith 
362273d9f13SBarry Smith   PetscFunctionBegin;
363273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
364273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g)
365273d9f13SBarry Smith   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
366273d9f13SBarry Smith   for (ii=1; ii<B->m; ii++) {
367273d9f13SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) {
368273d9f13SBarry Smith       SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]);
369273d9f13SBarry Smith     }
370273d9f13SBarry Smith   }
371273d9f13SBarry Smith   for (ii=0; ii<i[B->m]; ii++) {
372273d9f13SBarry Smith     if (j[ii] < 0 || j[ii] >= B->N) {
373273d9f13SBarry Smith       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
374273d9f13SBarry Smith     }
375273d9f13SBarry Smith   }
376273d9f13SBarry Smith #endif
377273d9f13SBarry Smith 
378273d9f13SBarry Smith   b->j      = j;
379273d9f13SBarry Smith   b->i      = i;
380273d9f13SBarry Smith   b->values = values;
381273d9f13SBarry Smith 
382273d9f13SBarry Smith   b->nz               = i[B->m];
383273d9f13SBarry Smith   b->diag             = 0;
384273d9f13SBarry Smith   b->symmetric        = PETSC_FALSE;
385273d9f13SBarry Smith   b->freeaij          = PETSC_TRUE;
386273d9f13SBarry Smith 
387273d9f13SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
388273d9f13SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
389273d9f13SBarry Smith   PetscFunctionReturn(0);
390273d9f13SBarry Smith }
391273d9f13SBarry Smith 
392b97cf49bSBarry Smith #undef __FUNC__
393b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIAdj"
394b97cf49bSBarry Smith /*@C
3953eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
396b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
397b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
398b97cf49bSBarry Smith 
399ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
400ef5ee4d1SLois Curfman McInnes 
401b97cf49bSBarry Smith    Input Parameters:
402c2e958feSBarry Smith +  comm - MPI communicator
4030752156aSBarry Smith .  m - number of local rows
404b97cf49bSBarry Smith .  n - number of columns
405b97cf49bSBarry Smith .  i - the indices into j for the start of each row
406329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
407ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
408329f5518SBarry Smith -  values -[optional] edge weights
409b97cf49bSBarry Smith 
410b97cf49bSBarry Smith    Output Parameter:
411b97cf49bSBarry Smith .  A - the matrix
412b97cf49bSBarry Smith 
41315091d37SBarry Smith    Level: intermediate
41415091d37SBarry Smith 
4154bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4164bc6d8c0SBarry Smith    MatSetValues().
417c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
4181198db28SBarry Smith    when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you
4191198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
420327e1a73SBarry Smith    Should not include the matrix diagonals.
421b97cf49bSBarry Smith 
422ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
423b97cf49bSBarry Smith 
42491e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering()
425b97cf49bSBarry Smith @*/
4263eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
427b97cf49bSBarry Smith {
428273d9f13SBarry Smith   int        ierr;
429b97cf49bSBarry Smith 
430433994e6SBarry Smith   PetscFunctionBegin;
431273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
432273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
433273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
4343a40ed3dSBarry Smith   PetscFunctionReturn(0);
435b97cf49bSBarry Smith }
436b97cf49bSBarry Smith 
437273d9f13SBarry Smith EXTERN_C_BEGIN
438c2e958feSBarry Smith #undef __FUNC__
439273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatConvertTo_MPIAdj"
440273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B)
441c2e958feSBarry Smith {
4424c49b128SBarry Smith   int      i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
443c2e958feSBarry Smith   Scalar   *ra;
444c2e958feSBarry Smith   MPI_Comm comm;
445c2e958feSBarry Smith 
446c2e958feSBarry Smith   PetscFunctionBegin;
447c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
448c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
449c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
450c2e958feSBarry Smith 
451c2e958feSBarry Smith   /* count the number of nonzeros per row */
452c2e958feSBarry Smith   for (i=0; i<m; i++) {
453c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
454c2e958feSBarry Smith     for (j=0; j<len; j++) {
455c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
456c2e958feSBarry Smith     }
457c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
458c2e958feSBarry Smith     nzeros += len;
459c2e958feSBarry Smith   }
460c2e958feSBarry Smith 
461c2e958feSBarry Smith   /* malloc space for nonzeros */
462c2e958feSBarry Smith   a  = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(a);
463c2e958feSBarry Smith   ia = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(ia);
464c2e958feSBarry Smith   ja = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(ja);
465c2e958feSBarry Smith 
466c2e958feSBarry Smith   nzeros = 0;
467c2e958feSBarry Smith   ia[0]  = 0;
468c2e958feSBarry Smith   for (i=0; i<m; i++) {
469c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
470c2e958feSBarry Smith     cnt     = 0;
471c2e958feSBarry Smith     for (j=0; j<len; j++) {
472c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
4734c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
474c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
475c2e958feSBarry Smith       }
476c2e958feSBarry Smith     }
477c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
478c2e958feSBarry Smith     nzeros += cnt;
479c2e958feSBarry Smith     ia[i+1] = nzeros;
480c2e958feSBarry Smith   }
481c2e958feSBarry Smith 
482c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
4833eda8832SBarry Smith   ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr);
484c2e958feSBarry Smith 
485c2e958feSBarry Smith   PetscFunctionReturn(0);
486c2e958feSBarry Smith }
487273d9f13SBarry Smith EXTERN_C_END
488433994e6SBarry Smith 
489433994e6SBarry Smith 
490433994e6SBarry Smith 
491433994e6SBarry Smith 
492433994e6SBarry Smith 
493