xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 0b3b1487d4e2d965bae3bb8a768736965a2af5d0)
1b97cf49bSBarry Smith #ifdef PETSC_RCS_HEADER
2*0b3b1487SBarry Smith static char vcid[] = "$Id: mpiadj.c,v 1.14 1998/06/18 15:15:56 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   }
74b97cf49bSBarry Smith 
753a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
7694d884c6SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
77b97cf49bSBarry Smith #endif
78b97cf49bSBarry Smith   if (a->diag) PetscFree(a->diag);
79b97cf49bSBarry Smith   PetscFree(a->i);
80b97cf49bSBarry Smith   PetscFree(a->j);
810752156aSBarry Smith   PetscFree(a->rowners);
82b97cf49bSBarry Smith   PetscFree(a);
83b97cf49bSBarry Smith 
8494d884c6SBarry Smith   PLogObjectDestroy(mat);
8594d884c6SBarry Smith   PetscHeaderDestroy(mat);
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
87b97cf49bSBarry Smith }
88b97cf49bSBarry Smith 
89b97cf49bSBarry Smith 
90b97cf49bSBarry Smith #undef __FUNC__
910752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj"
920752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
93b97cf49bSBarry Smith {
940752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
95b97cf49bSBarry Smith 
96b97cf49bSBarry Smith   if (op == MAT_STRUCTURALLY_SYMMETRIC) {
97b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
98b97cf49bSBarry Smith   } else {
99981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
100b97cf49bSBarry Smith   }
1013a40ed3dSBarry Smith   PetscFunctionReturn(0);
102b97cf49bSBarry Smith }
103b97cf49bSBarry Smith 
104b97cf49bSBarry Smith 
105b97cf49bSBarry Smith /*
106b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
107b97cf49bSBarry Smith */
108b97cf49bSBarry Smith 
109b97cf49bSBarry Smith #undef __FUNC__
1100752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj"
1110752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A)
112b97cf49bSBarry Smith {
1130752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
114b97cf49bSBarry Smith   int        i,j, *diag, m = a->m;
115b97cf49bSBarry Smith 
116b97cf49bSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
117b97cf49bSBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
118b97cf49bSBarry Smith   for ( i=0; i<a->m; i++ ) {
119b97cf49bSBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
120b97cf49bSBarry Smith       if (a->j[j] == i) {
121b97cf49bSBarry Smith         diag[i] = j;
122b97cf49bSBarry Smith         break;
123b97cf49bSBarry Smith       }
124b97cf49bSBarry Smith     }
125b97cf49bSBarry Smith   }
126b97cf49bSBarry Smith   a->diag = diag;
1273a40ed3dSBarry Smith   PetscFunctionReturn(0);
128b97cf49bSBarry Smith }
129b97cf49bSBarry Smith 
130b97cf49bSBarry Smith #undef __FUNC__
1310752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1320752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n)
133b97cf49bSBarry Smith {
1340752156aSBarry Smith   if (m) *m = A->M;
1350752156aSBarry Smith   if (n) *n = A->N;
1363a40ed3dSBarry Smith   PetscFunctionReturn(0);
137b97cf49bSBarry Smith }
138b97cf49bSBarry Smith 
139b97cf49bSBarry Smith #undef __FUNC__
1400752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1410752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n)
142b97cf49bSBarry Smith {
1430752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1440752156aSBarry Smith   if (m) *m = a->m;
1450752156aSBarry Smith   if (n) *n = A->N;
1463a40ed3dSBarry Smith   PetscFunctionReturn(0);
147b97cf49bSBarry Smith }
148b97cf49bSBarry Smith 
149b97cf49bSBarry Smith #undef __FUNC__
1500752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj"
1510752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
152b97cf49bSBarry Smith {
1530752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1540752156aSBarry Smith   *m = a->rstart; *n = a->rend;
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
156b97cf49bSBarry Smith }
157b97cf49bSBarry Smith 
158b97cf49bSBarry Smith #undef __FUNC__
1590752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj"
1600752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
161b97cf49bSBarry Smith {
1620752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
163b97cf49bSBarry Smith   int        *itmp;
164b97cf49bSBarry Smith 
1650752156aSBarry Smith   row -= a->rstart;
1660752156aSBarry Smith 
167a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
168b97cf49bSBarry Smith 
169b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
170b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
171b97cf49bSBarry Smith   if (idx) {
172b97cf49bSBarry Smith     itmp = a->j + a->i[row];
173b97cf49bSBarry Smith     if (*nz) {
174b97cf49bSBarry Smith       *idx = itmp;
175b97cf49bSBarry Smith     }
176b97cf49bSBarry Smith     else *idx = 0;
177b97cf49bSBarry Smith   }
1783a40ed3dSBarry Smith   PetscFunctionReturn(0);
179b97cf49bSBarry Smith }
180b97cf49bSBarry Smith 
181b97cf49bSBarry Smith #undef __FUNC__
1820752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj"
1830752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
184b97cf49bSBarry Smith {
1853a40ed3dSBarry Smith   PetscFunctionReturn(0);
186b97cf49bSBarry Smith }
187b97cf49bSBarry Smith 
188b97cf49bSBarry Smith #undef __FUNC__
1890752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj"
1900752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs)
191b97cf49bSBarry Smith {
192b97cf49bSBarry Smith   *bs = 1;
1933a40ed3dSBarry Smith   PetscFunctionReturn(0);
194b97cf49bSBarry Smith }
195b97cf49bSBarry Smith 
196b97cf49bSBarry Smith 
197b97cf49bSBarry Smith #undef __FUNC__
1980752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj"
1990752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg)
200b97cf49bSBarry Smith {
2010752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
202ca161407SBarry Smith  int         flag = 1,ierr;
203b97cf49bSBarry Smith 
204a8c6a408SBarry Smith   if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
205b97cf49bSBarry Smith 
206b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal, or no of nonzeros */
2070752156aSBarry Smith   if ((a->m != b->m ) ||( a->nz != b->nz)) {
2080752156aSBarry Smith     flag = 0;
209b97cf49bSBarry Smith   }
210b97cf49bSBarry Smith 
211b97cf49bSBarry Smith   /* if the a->i are the same */
212b97cf49bSBarry Smith   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
2130752156aSBarry Smith     flag = 0;
214b97cf49bSBarry Smith   }
215b97cf49bSBarry Smith 
216b97cf49bSBarry Smith   /* if a->j are the same */
217b97cf49bSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
2180752156aSBarry Smith     flag = 0;
219b97cf49bSBarry Smith   }
220b97cf49bSBarry Smith 
221ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
2220752156aSBarry Smith 
2230752156aSBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionReturn(0);
225b97cf49bSBarry Smith }
226b97cf49bSBarry Smith 
227b97cf49bSBarry Smith 
228b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
229*0b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
230*0b3b1487SBarry Smith        MatGetRow_MPIAdj,
231*0b3b1487SBarry Smith        MatRestoreRow_MPIAdj,
232b97cf49bSBarry Smith        0,
233b97cf49bSBarry Smith        0,
234b97cf49bSBarry Smith        0,
235b97cf49bSBarry Smith        0,
236*0b3b1487SBarry Smith        0,
237*0b3b1487SBarry Smith        0,
238*0b3b1487SBarry Smith        0,
239*0b3b1487SBarry Smith        0,
240*0b3b1487SBarry Smith        0,
241*0b3b1487SBarry Smith        0,
242*0b3b1487SBarry Smith        0,
243*0b3b1487SBarry Smith        0,
244*0b3b1487SBarry Smith        0,
245*0b3b1487SBarry Smith        MatEqual_MPIAdj,
246*0b3b1487SBarry Smith        0,
247*0b3b1487SBarry Smith        0,
248*0b3b1487SBarry Smith        0,
249*0b3b1487SBarry Smith        0,
250*0b3b1487SBarry Smith        0,
251*0b3b1487SBarry Smith        0,
252*0b3b1487SBarry Smith        MatSetOption_MPIAdj,
253*0b3b1487SBarry Smith        0,
254*0b3b1487SBarry Smith        0,
255*0b3b1487SBarry Smith        0,
256*0b3b1487SBarry Smith        0,
257*0b3b1487SBarry Smith        0,
258*0b3b1487SBarry Smith        0,
259*0b3b1487SBarry Smith        MatGetSize_MPIAdj,
260*0b3b1487SBarry Smith        MatGetLocalSize_MPIAdj,
261*0b3b1487SBarry Smith        MatGetOwnershipRange_MPIAdj,
262*0b3b1487SBarry Smith        0,
263*0b3b1487SBarry Smith        0,
264*0b3b1487SBarry Smith        0,
265*0b3b1487SBarry Smith        0,
266*0b3b1487SBarry Smith        0,
267*0b3b1487SBarry Smith        0,
268*0b3b1487SBarry Smith        0,
269*0b3b1487SBarry Smith        0,
270*0b3b1487SBarry Smith        0,
271*0b3b1487SBarry Smith        0,
272*0b3b1487SBarry Smith        0,
273*0b3b1487SBarry Smith        0,
274*0b3b1487SBarry Smith        0,
275*0b3b1487SBarry Smith        0,
276*0b3b1487SBarry Smith        0,
277*0b3b1487SBarry Smith        0,
278*0b3b1487SBarry Smith        0,
279*0b3b1487SBarry Smith        0,
280b97cf49bSBarry Smith        0,
2810752156aSBarry Smith        MatGetBlockSize_MPIAdj,
282b97cf49bSBarry Smith        0,
283b97cf49bSBarry Smith        0,
284b97cf49bSBarry Smith        0,
285b97cf49bSBarry Smith        0,
286b97cf49bSBarry Smith        0,
2870752156aSBarry Smith        0,
2880752156aSBarry Smith        0,
289*0b3b1487SBarry Smith        0,
290*0b3b1487SBarry Smith        0,
291*0b3b1487SBarry Smith        0,
292*0b3b1487SBarry Smith        0,
293*0b3b1487SBarry Smith        0,
294*0b3b1487SBarry Smith        MatGetMaps_Petsc};
295b97cf49bSBarry Smith 
296b97cf49bSBarry Smith 
297b97cf49bSBarry Smith #undef __FUNC__
2980752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj"
299b97cf49bSBarry Smith /*@C
3000752156aSBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
301b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
302b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
303b97cf49bSBarry Smith 
304ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
305ef5ee4d1SLois Curfman McInnes 
306b97cf49bSBarry Smith    Input Parameters:
307ef5ee4d1SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
3080752156aSBarry Smith .  m - number of local rows
309b97cf49bSBarry Smith .  n - number of columns
310b97cf49bSBarry Smith .  i - the indices into j for the start of each row
311ef5ee4d1SLois Curfman McInnes -  j - the column indices for each row (sorted for each row).
312ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
313b97cf49bSBarry Smith 
314b97cf49bSBarry Smith    Output Parameter:
315b97cf49bSBarry Smith .  A - the matrix
316b97cf49bSBarry Smith 
3174bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
3184bc6d8c0SBarry Smith    MatSetValues().
3194bc6d8c0SBarry Smith    You must NOT free the ii and jj arrays yourself. PETSc will free them
320b97cf49bSBarry Smith    when the matrix is destroyed.
321b97cf49bSBarry Smith 
322ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
323b97cf49bSBarry Smith 
3240752156aSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetReordering()
325b97cf49bSBarry Smith @*/
3260752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A)
327b97cf49bSBarry Smith {
328b97cf49bSBarry Smith   Mat        B;
3290752156aSBarry Smith   Mat_MPIAdj *b;
3300752156aSBarry Smith   int        ii,ierr, flg,size,rank;
331b97cf49bSBarry Smith 
332b97cf49bSBarry Smith   MPI_Comm_size(comm,&size);
3330752156aSBarry Smith   MPI_Comm_rank(comm,&rank);
334b97cf49bSBarry Smith 
335b97cf49bSBarry Smith   *A                  = 0;
336f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,comm,MatDestroy,MatView);
337b97cf49bSBarry Smith   PLogObjectCreate(B);
3380752156aSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_MPIAdj)); CHKPTRQ(b);
3390752156aSBarry Smith   PetscMemzero(b,sizeof(Mat_MPIAdj));
340*0b3b1487SBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
341e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_MPIAdj;
342e1311b90SBarry Smith   B->ops->view             = MatView_MPIAdj;
343b97cf49bSBarry Smith   B->factor           = 0;
344b97cf49bSBarry Smith   B->lupivotthreshold = 1.0;
345b97cf49bSBarry Smith   B->mapping          = 0;
346b97cf49bSBarry Smith   B->assembled        = PETSC_FALSE;
347b97cf49bSBarry Smith 
3480752156aSBarry Smith   b->m = m; B->m = m;
349ca161407SBarry Smith   ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
3500752156aSBarry Smith   B->n = n; B->N = n;
3510752156aSBarry Smith 
3520752156aSBarry Smith   b->rowners = (int *) PetscMalloc((size+1)*sizeof(int)); CHKPTRQ(b->rowners);
3530752156aSBarry Smith   PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
354ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
3550752156aSBarry Smith   b->rowners[0] = 0;
3560752156aSBarry Smith   for ( ii=2; ii<=size; ii++ ) {
3570752156aSBarry Smith     b->rowners[ii] += b->rowners[ii-1];
3580752156aSBarry Smith   }
3590752156aSBarry Smith   b->rstart = b->rowners[rank];
3600752156aSBarry Smith   b->rend   = b->rowners[rank+1];
361b97cf49bSBarry Smith 
362b97cf49bSBarry Smith   b->j  = j;
363b97cf49bSBarry Smith   b->i  = i;
364b97cf49bSBarry Smith 
365b97cf49bSBarry Smith   b->nz               = i[m];
366b97cf49bSBarry Smith   b->diag             = 0;
367b97cf49bSBarry Smith   b->symmetric        = PETSC_FALSE;
368b97cf49bSBarry Smith 
369b97cf49bSBarry Smith   *A = B;
370b97cf49bSBarry Smith 
371b97cf49bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
372b97cf49bSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
373b97cf49bSBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
374b97cf49bSBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3753a40ed3dSBarry Smith   PetscFunctionReturn(0);
376b97cf49bSBarry Smith }
377b97cf49bSBarry Smith 
378b97cf49bSBarry Smith 
379b97cf49bSBarry Smith 
380