xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 488ecbaffa467bb032d31d7eb20bc6d0ef6d986c)
1b97cf49bSBarry Smith #ifdef PETSC_RCS_HEADER
2*488ecbafSBarry Smith static char vcid[] = "$Id: mpiadj.c,v 1.16 1998/07/14 03:16:14 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 
9b97cf49bSBarry Smith #include "pinclude/pviewer.h"
10b97cf49bSBarry Smith #include "sys.h"
110752156aSBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
12b97cf49bSBarry Smith 
13b97cf49bSBarry Smith 
14b97cf49bSBarry Smith #undef __FUNC__
150752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj_ASCII"
160752156aSBarry Smith extern int MatView_MPIAdj_ASCII(Mat A,Viewer viewer)
17b97cf49bSBarry Smith {
180752156aSBarry Smith   Mat_MPIAdj  *a = (Mat_MPIAdj *) A->data;
19b97cf49bSBarry Smith   int         ierr, i,j, m = a->m,  format;
20b97cf49bSBarry Smith   FILE        *fd;
21b97cf49bSBarry Smith   char        *outputname;
220752156aSBarry Smith   MPI_Comm    comm = A->comm;
23b97cf49bSBarry Smith 
24b97cf49bSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
25b97cf49bSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
26b97cf49bSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
27b97cf49bSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
283a40ed3dSBarry Smith     PetscFunctionReturn(0);
290752156aSBarry Smith   } else {
30b97cf49bSBarry Smith     for ( i=0; i<m; i++ ) {
310752156aSBarry Smith       PetscSynchronizedFPrintf(comm,fd,"row %d:",i+a->rstart);
32b97cf49bSBarry Smith       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
330752156aSBarry Smith         PetscSynchronizedFPrintf(comm,fd," %d ",a->j[j]);
340752156aSBarry Smith       }
350752156aSBarry Smith       PetscSynchronizedFPrintf(comm,fd,"\n");
36b97cf49bSBarry Smith     }
37b97cf49bSBarry Smith   }
380752156aSBarry Smith   PetscSynchronizedFlush(comm);
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
40b97cf49bSBarry Smith }
41b97cf49bSBarry Smith 
42b97cf49bSBarry Smith #undef __FUNC__
430752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj"
44e1311b90SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer)
45b97cf49bSBarry Smith {
46b97cf49bSBarry Smith   ViewerType  vtype;
47b97cf49bSBarry Smith   int         ierr;
48b97cf49bSBarry Smith 
49b97cf49bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
500752156aSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
513a40ed3dSBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
525cd90555SBarry Smith   } else {
535cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
54b97cf49bSBarry Smith   }
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56b97cf49bSBarry Smith }
57b97cf49bSBarry Smith 
58b97cf49bSBarry Smith #undef __FUNC__
590752156aSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj"
6094d884c6SBarry Smith int MatDestroy_MPIAdj(Mat mat)
61b97cf49bSBarry Smith {
6294d884c6SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) mat->data;
6394d884c6SBarry Smith   int        ierr;
6494d884c6SBarry Smith 
6594d884c6SBarry Smith   PetscFunctionBegin;
6694d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
6794d884c6SBarry Smith 
6894d884c6SBarry Smith   if (mat->mapping) {
6994d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
7094d884c6SBarry Smith   }
7194d884c6SBarry Smith   if (mat->bmapping) {
7294d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
7394d884c6SBarry Smith   }
74c7fcc2eaSBarry Smith   if (mat->rmap) {
75c7fcc2eaSBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
76c7fcc2eaSBarry Smith   }
77c7fcc2eaSBarry Smith   if (mat->cmap) {
78c7fcc2eaSBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
79c7fcc2eaSBarry Smith   }
80b97cf49bSBarry Smith 
813a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
8294d884c6SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
83b97cf49bSBarry Smith #endif
84b97cf49bSBarry Smith   if (a->diag) PetscFree(a->diag);
85b97cf49bSBarry Smith   PetscFree(a->i);
86b97cf49bSBarry Smith   PetscFree(a->j);
870752156aSBarry Smith   PetscFree(a->rowners);
88b97cf49bSBarry Smith   PetscFree(a);
89b97cf49bSBarry Smith 
9094d884c6SBarry Smith   PLogObjectDestroy(mat);
9194d884c6SBarry Smith   PetscHeaderDestroy(mat);
923a40ed3dSBarry Smith   PetscFunctionReturn(0);
93b97cf49bSBarry Smith }
94b97cf49bSBarry Smith 
95b97cf49bSBarry Smith 
96b97cf49bSBarry Smith #undef __FUNC__
970752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj"
980752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
99b97cf49bSBarry Smith {
1000752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
101b97cf49bSBarry Smith 
102b97cf49bSBarry Smith   if (op == MAT_STRUCTURALLY_SYMMETRIC) {
103b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
104b97cf49bSBarry Smith   } else {
105981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
106b97cf49bSBarry Smith   }
1073a40ed3dSBarry Smith   PetscFunctionReturn(0);
108b97cf49bSBarry Smith }
109b97cf49bSBarry Smith 
110b97cf49bSBarry Smith 
111b97cf49bSBarry Smith /*
112b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
113b97cf49bSBarry Smith */
114b97cf49bSBarry Smith 
115b97cf49bSBarry Smith #undef __FUNC__
1160752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj"
1170752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A)
118b97cf49bSBarry Smith {
1190752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
120b97cf49bSBarry Smith   int        i,j, *diag, m = a->m;
121b97cf49bSBarry Smith 
122b97cf49bSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
123b97cf49bSBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
124b97cf49bSBarry Smith   for ( i=0; i<a->m; i++ ) {
125b97cf49bSBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
126b97cf49bSBarry Smith       if (a->j[j] == i) {
127b97cf49bSBarry Smith         diag[i] = j;
128b97cf49bSBarry Smith         break;
129b97cf49bSBarry Smith       }
130b97cf49bSBarry Smith     }
131b97cf49bSBarry Smith   }
132b97cf49bSBarry Smith   a->diag = diag;
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
134b97cf49bSBarry Smith }
135b97cf49bSBarry Smith 
136b97cf49bSBarry Smith #undef __FUNC__
1370752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1380752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n)
139b97cf49bSBarry Smith {
1400752156aSBarry Smith   if (m) *m = A->M;
1410752156aSBarry Smith   if (n) *n = A->N;
1423a40ed3dSBarry Smith   PetscFunctionReturn(0);
143b97cf49bSBarry Smith }
144b97cf49bSBarry Smith 
145b97cf49bSBarry Smith #undef __FUNC__
1460752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1470752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n)
148b97cf49bSBarry Smith {
1490752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1500752156aSBarry Smith   if (m) *m = a->m;
1510752156aSBarry Smith   if (n) *n = A->N;
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
153b97cf49bSBarry Smith }
154b97cf49bSBarry Smith 
155b97cf49bSBarry Smith #undef __FUNC__
1560752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj"
1570752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
158b97cf49bSBarry Smith {
1590752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1600752156aSBarry Smith   *m = a->rstart; *n = a->rend;
1613a40ed3dSBarry Smith   PetscFunctionReturn(0);
162b97cf49bSBarry Smith }
163b97cf49bSBarry Smith 
164b97cf49bSBarry Smith #undef __FUNC__
1650752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj"
1660752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
167b97cf49bSBarry Smith {
1680752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
169b97cf49bSBarry Smith   int        *itmp;
170b97cf49bSBarry Smith 
1710752156aSBarry Smith   row -= a->rstart;
1720752156aSBarry Smith 
173a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
174b97cf49bSBarry Smith 
175b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
176b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
177b97cf49bSBarry Smith   if (idx) {
178b97cf49bSBarry Smith     itmp = a->j + a->i[row];
179b97cf49bSBarry Smith     if (*nz) {
180b97cf49bSBarry Smith       *idx = itmp;
181b97cf49bSBarry Smith     }
182b97cf49bSBarry Smith     else *idx = 0;
183b97cf49bSBarry Smith   }
1843a40ed3dSBarry Smith   PetscFunctionReturn(0);
185b97cf49bSBarry Smith }
186b97cf49bSBarry Smith 
187b97cf49bSBarry Smith #undef __FUNC__
1880752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj"
1890752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
190b97cf49bSBarry Smith {
1913a40ed3dSBarry Smith   PetscFunctionReturn(0);
192b97cf49bSBarry Smith }
193b97cf49bSBarry Smith 
194b97cf49bSBarry Smith #undef __FUNC__
1950752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj"
1960752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs)
197b97cf49bSBarry Smith {
198b97cf49bSBarry Smith   *bs = 1;
1993a40ed3dSBarry Smith   PetscFunctionReturn(0);
200b97cf49bSBarry Smith }
201b97cf49bSBarry Smith 
202b97cf49bSBarry Smith 
203b97cf49bSBarry Smith #undef __FUNC__
2040752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj"
2050752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg)
206b97cf49bSBarry Smith {
2070752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
208ca161407SBarry Smith  int         flag = 1,ierr;
209b97cf49bSBarry Smith 
210a8c6a408SBarry Smith   if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
211b97cf49bSBarry Smith 
212b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal, or no of nonzeros */
2130752156aSBarry Smith   if ((a->m != b->m ) ||( a->nz != b->nz)) {
2140752156aSBarry Smith     flag = 0;
215b97cf49bSBarry Smith   }
216b97cf49bSBarry Smith 
217b97cf49bSBarry Smith   /* if the a->i are the same */
218b97cf49bSBarry Smith   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
2190752156aSBarry Smith     flag = 0;
220b97cf49bSBarry Smith   }
221b97cf49bSBarry Smith 
222b97cf49bSBarry Smith   /* if a->j are the same */
223b97cf49bSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
2240752156aSBarry Smith     flag = 0;
225b97cf49bSBarry Smith   }
226b97cf49bSBarry Smith 
227ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
2280752156aSBarry Smith 
2290752156aSBarry Smith 
2303a40ed3dSBarry Smith   PetscFunctionReturn(0);
231b97cf49bSBarry Smith }
232b97cf49bSBarry Smith 
233b97cf49bSBarry Smith 
234b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2350b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2360b3b1487SBarry Smith        MatGetRow_MPIAdj,
2370b3b1487SBarry Smith        MatRestoreRow_MPIAdj,
238b97cf49bSBarry Smith        0,
239b97cf49bSBarry Smith        0,
240b97cf49bSBarry Smith        0,
241b97cf49bSBarry Smith        0,
2420b3b1487SBarry Smith        0,
2430b3b1487SBarry 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        MatEqual_MPIAdj,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        MatSetOption_MPIAdj,
2590b3b1487SBarry Smith        0,
2600b3b1487SBarry Smith        0,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
2640b3b1487SBarry Smith        0,
2650b3b1487SBarry Smith        MatGetSize_MPIAdj,
2660b3b1487SBarry Smith        MatGetLocalSize_MPIAdj,
2670b3b1487SBarry Smith        MatGetOwnershipRange_MPIAdj,
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,
2830b3b1487SBarry Smith        0,
2840b3b1487SBarry Smith        0,
2850b3b1487SBarry Smith        0,
286b97cf49bSBarry Smith        0,
2870752156aSBarry Smith        MatGetBlockSize_MPIAdj,
288b97cf49bSBarry Smith        0,
289b97cf49bSBarry Smith        0,
290b97cf49bSBarry Smith        0,
291b97cf49bSBarry Smith        0,
292b97cf49bSBarry Smith        0,
2930752156aSBarry Smith        0,
2940752156aSBarry Smith        0,
2950b3b1487SBarry Smith        0,
2960b3b1487SBarry Smith        0,
2970b3b1487SBarry Smith        0,
2980b3b1487SBarry Smith        0,
2990b3b1487SBarry Smith        0,
3000b3b1487SBarry Smith        MatGetMaps_Petsc};
301b97cf49bSBarry Smith 
302b97cf49bSBarry Smith 
303b97cf49bSBarry Smith #undef __FUNC__
3040752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj"
305b97cf49bSBarry Smith /*@C
3060752156aSBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
307b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
308b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
309b97cf49bSBarry Smith 
310ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
311ef5ee4d1SLois Curfman McInnes 
312b97cf49bSBarry Smith    Input Parameters:
313ef5ee4d1SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
3140752156aSBarry Smith .  m - number of local rows
315b97cf49bSBarry Smith .  n - number of columns
316b97cf49bSBarry Smith .  i - the indices into j for the start of each row
317ef5ee4d1SLois Curfman McInnes -  j - the column indices for each row (sorted for each row).
318ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
319b97cf49bSBarry Smith 
320b97cf49bSBarry Smith    Output Parameter:
321b97cf49bSBarry Smith .  A - the matrix
322b97cf49bSBarry Smith 
3234bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
3244bc6d8c0SBarry Smith    MatSetValues().
3254bc6d8c0SBarry Smith    You must NOT free the ii and jj arrays yourself. PETSc will free them
326b97cf49bSBarry Smith    when the matrix is destroyed.
327b97cf49bSBarry Smith 
328ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
329b97cf49bSBarry Smith 
3300752156aSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetReordering()
331b97cf49bSBarry Smith @*/
3320752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A)
333b97cf49bSBarry Smith {
334b97cf49bSBarry Smith   Mat        B;
3350752156aSBarry Smith   Mat_MPIAdj *b;
3360752156aSBarry Smith   int        ii,ierr, flg,size,rank;
337b97cf49bSBarry Smith 
338b97cf49bSBarry Smith   MPI_Comm_size(comm,&size);
3390752156aSBarry Smith   MPI_Comm_rank(comm,&rank);
340b97cf49bSBarry Smith 
341b97cf49bSBarry Smith   *A                  = 0;
342f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,comm,MatDestroy,MatView);
343b97cf49bSBarry Smith   PLogObjectCreate(B);
3440752156aSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_MPIAdj)); CHKPTRQ(b);
3450752156aSBarry Smith   PetscMemzero(b,sizeof(Mat_MPIAdj));
3460b3b1487SBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
347e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_MPIAdj;
348e1311b90SBarry Smith   B->ops->view             = MatView_MPIAdj;
349b97cf49bSBarry Smith   B->factor           = 0;
350b97cf49bSBarry Smith   B->lupivotthreshold = 1.0;
351b97cf49bSBarry Smith   B->mapping          = 0;
352b97cf49bSBarry Smith   B->assembled        = PETSC_FALSE;
353b97cf49bSBarry Smith 
3540752156aSBarry Smith   b->m = m; B->m = m;
355ca161407SBarry Smith   ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
3560752156aSBarry Smith   B->n = n; B->N = n;
3570752156aSBarry Smith 
358c7fcc2eaSBarry Smith   /* the information in the maps duplicates the information computed below, eventually
359c7fcc2eaSBarry Smith      we should remove the duplicate information that is not contained in the maps */
360*488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,B->M,&B->rmap);CHKERRQ(ierr);
361c7fcc2eaSBarry Smith   /* we don't know the "local columns" so just use the row information :-( */
362*488ecbafSBarry Smith   ierr = MapCreateMPI(comm,m,B->M,&B->cmap);CHKERRQ(ierr);
363c7fcc2eaSBarry Smith 
3640752156aSBarry Smith   b->rowners = (int *) PetscMalloc((size+1)*sizeof(int)); CHKPTRQ(b->rowners);
3650752156aSBarry Smith   PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
366ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
3670752156aSBarry Smith   b->rowners[0] = 0;
3680752156aSBarry Smith   for ( ii=2; ii<=size; ii++ ) {
3690752156aSBarry Smith     b->rowners[ii] += b->rowners[ii-1];
3700752156aSBarry Smith   }
3710752156aSBarry Smith   b->rstart = b->rowners[rank];
3720752156aSBarry Smith   b->rend   = b->rowners[rank+1];
373b97cf49bSBarry Smith 
374b97cf49bSBarry Smith   b->j  = j;
375b97cf49bSBarry Smith   b->i  = i;
376b97cf49bSBarry Smith 
377b97cf49bSBarry Smith   b->nz               = i[m];
378b97cf49bSBarry Smith   b->diag             = 0;
379b97cf49bSBarry Smith   b->symmetric        = PETSC_FALSE;
380b97cf49bSBarry Smith 
381b97cf49bSBarry Smith   *A = B;
382b97cf49bSBarry Smith 
383b97cf49bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
384b97cf49bSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
385b97cf49bSBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
386b97cf49bSBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
388b97cf49bSBarry Smith }
389b97cf49bSBarry Smith 
390b97cf49bSBarry Smith 
391b97cf49bSBarry Smith 
392