xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 3f1db9ec2fd39765c6c3a00831044586630c4cca)
1b97cf49bSBarry Smith #ifdef PETSC_RCS_HEADER
2*3f1db9ecSBarry Smith static char vcid[] = "$Id: mpiadj.c,v 1.18 1998/12/03 04:01:42 bsmith Exp bsmith $";
3b97cf49bSBarry Smith #endif
4b97cf49bSBarry Smith 
5b97cf49bSBarry Smith /*
6b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
7b97cf49bSBarry Smith */
8b97cf49bSBarry Smith #include "sys.h"
90752156aSBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
10b97cf49bSBarry Smith 
11b97cf49bSBarry Smith #undef __FUNC__
120752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj_ASCII"
130752156aSBarry Smith extern int MatView_MPIAdj_ASCII(Mat A,Viewer viewer)
14b97cf49bSBarry Smith {
150752156aSBarry Smith   Mat_MPIAdj  *a = (Mat_MPIAdj *) A->data;
16b97cf49bSBarry Smith   int         ierr, i,j, m = a->m,  format;
17b97cf49bSBarry Smith   FILE        *fd;
18b97cf49bSBarry Smith   char        *outputname;
190752156aSBarry Smith   MPI_Comm    comm = A->comm;
20b97cf49bSBarry Smith 
21b97cf49bSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
2277ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname); CHKERRQ(ierr);
23b97cf49bSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
24b97cf49bSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
253a40ed3dSBarry Smith     PetscFunctionReturn(0);
260752156aSBarry Smith   } else {
27b97cf49bSBarry Smith     for ( i=0; i<m; i++ ) {
280752156aSBarry Smith       PetscSynchronizedFPrintf(comm,fd,"row %d:",i+a->rstart);
29b97cf49bSBarry Smith       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
300752156aSBarry Smith         PetscSynchronizedFPrintf(comm,fd," %d ",a->j[j]);
310752156aSBarry Smith       }
320752156aSBarry Smith       PetscSynchronizedFPrintf(comm,fd,"\n");
33b97cf49bSBarry Smith     }
34b97cf49bSBarry Smith   }
350752156aSBarry Smith   PetscSynchronizedFlush(comm);
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
37b97cf49bSBarry Smith }
38b97cf49bSBarry Smith 
39b97cf49bSBarry Smith #undef __FUNC__
400752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj"
41e1311b90SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer)
42b97cf49bSBarry Smith {
43b97cf49bSBarry Smith   ViewerType  vtype;
44b97cf49bSBarry Smith   int         ierr;
45b97cf49bSBarry Smith 
46b97cf49bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
47*3f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)){
483a40ed3dSBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
495cd90555SBarry Smith   } else {
505cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
51b97cf49bSBarry Smith   }
523a40ed3dSBarry Smith   PetscFunctionReturn(0);
53b97cf49bSBarry Smith }
54b97cf49bSBarry Smith 
55b97cf49bSBarry Smith #undef __FUNC__
560752156aSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj"
5794d884c6SBarry Smith int MatDestroy_MPIAdj(Mat mat)
58b97cf49bSBarry Smith {
5994d884c6SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) mat->data;
6094d884c6SBarry Smith   int        ierr;
6194d884c6SBarry Smith 
6294d884c6SBarry Smith   PetscFunctionBegin;
6394d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
6494d884c6SBarry Smith 
6594d884c6SBarry Smith   if (mat->mapping) {
6694d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
6794d884c6SBarry Smith   }
6894d884c6SBarry Smith   if (mat->bmapping) {
6994d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
7094d884c6SBarry Smith   }
71c7fcc2eaSBarry Smith   if (mat->rmap) {
72c7fcc2eaSBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
73c7fcc2eaSBarry Smith   }
74c7fcc2eaSBarry Smith   if (mat->cmap) {
75c7fcc2eaSBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
76c7fcc2eaSBarry Smith   }
77b97cf49bSBarry Smith 
783a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
7994d884c6SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
80b97cf49bSBarry Smith #endif
81b97cf49bSBarry Smith   if (a->diag) PetscFree(a->diag);
82b97cf49bSBarry Smith   PetscFree(a->i);
83b97cf49bSBarry Smith   PetscFree(a->j);
840752156aSBarry Smith   PetscFree(a->rowners);
85b97cf49bSBarry Smith   PetscFree(a);
86b97cf49bSBarry Smith 
8794d884c6SBarry Smith   PLogObjectDestroy(mat);
8894d884c6SBarry Smith   PetscHeaderDestroy(mat);
893a40ed3dSBarry Smith   PetscFunctionReturn(0);
90b97cf49bSBarry Smith }
91b97cf49bSBarry Smith 
92b97cf49bSBarry Smith 
93b97cf49bSBarry Smith #undef __FUNC__
940752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj"
950752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
96b97cf49bSBarry Smith {
970752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
98b97cf49bSBarry Smith 
99b97cf49bSBarry Smith   if (op == MAT_STRUCTURALLY_SYMMETRIC) {
100b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
101b97cf49bSBarry Smith   } else {
102981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
103b97cf49bSBarry Smith   }
1043a40ed3dSBarry Smith   PetscFunctionReturn(0);
105b97cf49bSBarry Smith }
106b97cf49bSBarry Smith 
107b97cf49bSBarry Smith 
108b97cf49bSBarry Smith /*
109b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
110b97cf49bSBarry Smith */
111b97cf49bSBarry Smith 
112b97cf49bSBarry Smith #undef __FUNC__
1130752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj"
1140752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A)
115b97cf49bSBarry Smith {
1160752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
117b97cf49bSBarry Smith   int        i,j, *diag, m = a->m;
118b97cf49bSBarry Smith 
119b97cf49bSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
120b97cf49bSBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
121b97cf49bSBarry Smith   for ( i=0; i<a->m; i++ ) {
122b97cf49bSBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
123b97cf49bSBarry Smith       if (a->j[j] == i) {
124b97cf49bSBarry Smith         diag[i] = j;
125b97cf49bSBarry Smith         break;
126b97cf49bSBarry Smith       }
127b97cf49bSBarry Smith     }
128b97cf49bSBarry Smith   }
129b97cf49bSBarry Smith   a->diag = diag;
1303a40ed3dSBarry Smith   PetscFunctionReturn(0);
131b97cf49bSBarry Smith }
132b97cf49bSBarry Smith 
133b97cf49bSBarry Smith #undef __FUNC__
1340752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1350752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n)
136b97cf49bSBarry Smith {
1370752156aSBarry Smith   if (m) *m = A->M;
1380752156aSBarry Smith   if (n) *n = A->N;
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
140b97cf49bSBarry Smith }
141b97cf49bSBarry Smith 
142b97cf49bSBarry Smith #undef __FUNC__
1430752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1440752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n)
145b97cf49bSBarry Smith {
1460752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1470752156aSBarry Smith   if (m) *m = a->m;
1480752156aSBarry Smith   if (n) *n = A->N;
1493a40ed3dSBarry Smith   PetscFunctionReturn(0);
150b97cf49bSBarry Smith }
151b97cf49bSBarry Smith 
152b97cf49bSBarry Smith #undef __FUNC__
1530752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj"
1540752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
155b97cf49bSBarry Smith {
1560752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1570752156aSBarry Smith   *m = a->rstart; *n = a->rend;
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
159b97cf49bSBarry Smith }
160b97cf49bSBarry Smith 
161b97cf49bSBarry Smith #undef __FUNC__
1620752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj"
1630752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
164b97cf49bSBarry Smith {
1650752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
166b97cf49bSBarry Smith   int        *itmp;
167b97cf49bSBarry Smith 
1680752156aSBarry Smith   row -= a->rstart;
1690752156aSBarry Smith 
170a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
171b97cf49bSBarry Smith 
172b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
173b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
174b97cf49bSBarry Smith   if (idx) {
175b97cf49bSBarry Smith     itmp = a->j + a->i[row];
176b97cf49bSBarry Smith     if (*nz) {
177b97cf49bSBarry Smith       *idx = itmp;
178b97cf49bSBarry Smith     }
179b97cf49bSBarry Smith     else *idx = 0;
180b97cf49bSBarry Smith   }
1813a40ed3dSBarry Smith   PetscFunctionReturn(0);
182b97cf49bSBarry Smith }
183b97cf49bSBarry Smith 
184b97cf49bSBarry Smith #undef __FUNC__
1850752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj"
1860752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
187b97cf49bSBarry Smith {
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
189b97cf49bSBarry Smith }
190b97cf49bSBarry Smith 
191b97cf49bSBarry Smith #undef __FUNC__
1920752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj"
1930752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs)
194b97cf49bSBarry Smith {
195b97cf49bSBarry Smith   *bs = 1;
1963a40ed3dSBarry Smith   PetscFunctionReturn(0);
197b97cf49bSBarry Smith }
198b97cf49bSBarry Smith 
199b97cf49bSBarry Smith 
200b97cf49bSBarry Smith #undef __FUNC__
2010752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj"
2020752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg)
203b97cf49bSBarry Smith {
2040752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
205ca161407SBarry Smith  int         flag = 1,ierr;
206b97cf49bSBarry Smith 
207a8c6a408SBarry Smith   if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
208b97cf49bSBarry Smith 
209b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal, or no of nonzeros */
2100752156aSBarry Smith   if ((a->m != b->m ) ||( a->nz != b->nz)) {
2110752156aSBarry Smith     flag = 0;
212b97cf49bSBarry Smith   }
213b97cf49bSBarry Smith 
214b97cf49bSBarry Smith   /* if the a->i are the same */
215b97cf49bSBarry Smith   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
2160752156aSBarry Smith     flag = 0;
217b97cf49bSBarry Smith   }
218b97cf49bSBarry Smith 
219b97cf49bSBarry Smith   /* if a->j are the same */
220b97cf49bSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
2210752156aSBarry Smith     flag = 0;
222b97cf49bSBarry Smith   }
223b97cf49bSBarry Smith 
224ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
2250752156aSBarry Smith 
2260752156aSBarry Smith 
2273a40ed3dSBarry Smith   PetscFunctionReturn(0);
228b97cf49bSBarry Smith }
229b97cf49bSBarry Smith 
230b97cf49bSBarry Smith 
231b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2320b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2330b3b1487SBarry Smith        MatGetRow_MPIAdj,
2340b3b1487SBarry Smith        MatRestoreRow_MPIAdj,
235b97cf49bSBarry Smith        0,
236b97cf49bSBarry Smith        0,
237b97cf49bSBarry Smith        0,
238b97cf49bSBarry Smith        0,
2390b3b1487SBarry Smith        0,
2400b3b1487SBarry Smith        0,
2410b3b1487SBarry Smith        0,
2420b3b1487SBarry Smith        0,
2430b3b1487SBarry Smith        0,
2440b3b1487SBarry Smith        0,
2450b3b1487SBarry Smith        0,
2460b3b1487SBarry Smith        0,
2470b3b1487SBarry Smith        0,
2480b3b1487SBarry Smith        MatEqual_MPIAdj,
2490b3b1487SBarry Smith        0,
2500b3b1487SBarry Smith        0,
2510b3b1487SBarry Smith        0,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        MatSetOption_MPIAdj,
2560b3b1487SBarry Smith        0,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
2600b3b1487SBarry Smith        0,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        MatGetSize_MPIAdj,
2630b3b1487SBarry Smith        MatGetLocalSize_MPIAdj,
2640b3b1487SBarry Smith        MatGetOwnershipRange_MPIAdj,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
2670b3b1487SBarry Smith        0,
2680b3b1487SBarry 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,
283b97cf49bSBarry Smith        0,
2840752156aSBarry Smith        MatGetBlockSize_MPIAdj,
285b97cf49bSBarry Smith        0,
286b97cf49bSBarry Smith        0,
287b97cf49bSBarry Smith        0,
288b97cf49bSBarry Smith        0,
289b97cf49bSBarry Smith        0,
2900752156aSBarry Smith        0,
2910752156aSBarry Smith        0,
2920b3b1487SBarry Smith        0,
2930b3b1487SBarry Smith        0,
2940b3b1487SBarry Smith        0,
2950b3b1487SBarry Smith        0,
2960b3b1487SBarry Smith        0,
2970b3b1487SBarry Smith        MatGetMaps_Petsc};
298b97cf49bSBarry Smith 
299b97cf49bSBarry Smith 
300b97cf49bSBarry Smith #undef __FUNC__
3010752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj"
302b97cf49bSBarry Smith /*@C
3030752156aSBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
304b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
305b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
306b97cf49bSBarry Smith 
307ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
308ef5ee4d1SLois Curfman McInnes 
309b97cf49bSBarry Smith    Input Parameters:
310ef5ee4d1SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
3110752156aSBarry Smith .  m - number of local rows
312b97cf49bSBarry Smith .  n - number of columns
313b97cf49bSBarry Smith .  i - the indices into j for the start of each row
314ef5ee4d1SLois Curfman McInnes -  j - the column indices for each row (sorted for each row).
315ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
316b97cf49bSBarry Smith 
317b97cf49bSBarry Smith    Output Parameter:
318b97cf49bSBarry Smith .  A - the matrix
319b97cf49bSBarry Smith 
3204bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
3214bc6d8c0SBarry Smith    MatSetValues().
3224bc6d8c0SBarry Smith    You must NOT free the ii and jj arrays yourself. PETSc will free them
323b97cf49bSBarry Smith    when the matrix is destroyed.
324b97cf49bSBarry Smith 
325ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
326b97cf49bSBarry Smith 
3270752156aSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetReordering()
328b97cf49bSBarry Smith @*/
3290752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A)
330b97cf49bSBarry Smith {
331b97cf49bSBarry Smith   Mat        B;
3320752156aSBarry Smith   Mat_MPIAdj *b;
3330752156aSBarry Smith   int        ii,ierr, flg,size,rank;
334b97cf49bSBarry Smith 
335b97cf49bSBarry Smith   MPI_Comm_size(comm,&size);
3360752156aSBarry Smith   MPI_Comm_rank(comm,&rank);
337b97cf49bSBarry Smith 
338b97cf49bSBarry Smith   *A                  = 0;
339*3f1db9ecSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,"Mat",comm,MatDestroy,MatView);
340b97cf49bSBarry Smith   PLogObjectCreate(B);
3410752156aSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_MPIAdj)); CHKPTRQ(b);
3420752156aSBarry Smith   PetscMemzero(b,sizeof(Mat_MPIAdj));
3430b3b1487SBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
344e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_MPIAdj;
345e1311b90SBarry Smith   B->ops->view             = MatView_MPIAdj;
346b97cf49bSBarry Smith   B->factor           = 0;
347b97cf49bSBarry Smith   B->lupivotthreshold = 1.0;
348b97cf49bSBarry Smith   B->mapping          = 0;
349b97cf49bSBarry Smith   B->assembled        = PETSC_FALSE;
350b97cf49bSBarry Smith 
3510752156aSBarry Smith   b->m = m; B->m = m;
352ca161407SBarry Smith   ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
3530752156aSBarry Smith   B->n = n; B->N = n;
3540752156aSBarry Smith 
355c7fcc2eaSBarry Smith   /* the information in the maps duplicates the information computed below, eventually
356c7fcc2eaSBarry Smith      we should remove the duplicate information that is not contained in the maps */
357488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,B->M,&B->rmap);CHKERRQ(ierr);
358c7fcc2eaSBarry Smith   /* we don't know the "local columns" so just use the row information :-( */
359488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,B->M,&B->cmap);CHKERRQ(ierr);
360c7fcc2eaSBarry Smith 
3610752156aSBarry Smith   b->rowners = (int *) PetscMalloc((size+1)*sizeof(int)); CHKPTRQ(b->rowners);
3620752156aSBarry Smith   PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
363ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
3640752156aSBarry Smith   b->rowners[0] = 0;
3650752156aSBarry Smith   for ( ii=2; ii<=size; ii++ ) {
3660752156aSBarry Smith     b->rowners[ii] += b->rowners[ii-1];
3670752156aSBarry Smith   }
3680752156aSBarry Smith   b->rstart = b->rowners[rank];
3690752156aSBarry Smith   b->rend   = b->rowners[rank+1];
370b97cf49bSBarry Smith 
371b97cf49bSBarry Smith   b->j  = j;
372b97cf49bSBarry Smith   b->i  = i;
373b97cf49bSBarry Smith 
374b97cf49bSBarry Smith   b->nz               = i[m];
375b97cf49bSBarry Smith   b->diag             = 0;
376b97cf49bSBarry Smith   b->symmetric        = PETSC_FALSE;
377b97cf49bSBarry Smith 
378b97cf49bSBarry Smith   *A = B;
379b97cf49bSBarry Smith 
380b97cf49bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
381b97cf49bSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
382b97cf49bSBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
383b97cf49bSBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3843a40ed3dSBarry Smith   PetscFunctionReturn(0);
385b97cf49bSBarry Smith }
386b97cf49bSBarry Smith 
387b97cf49bSBarry Smith 
388b97cf49bSBarry Smith 
389