xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision bd019fc1a6f956260a0ee235b3d504ca5a0b1f2a)
159557b74SHong Zhang 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
4c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
559557b74SHong Zhang 
6cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
74e5e7fe4SHong Zhang {
84e5e7fe4SHong Zhang   Mat            B;
94e5e7fe4SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
104e5e7fe4SHong Zhang   Mat_SeqAIJ     *b;
11dfbe8321SBarry Smith   PetscErrorCode ierr;
12d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*rowlengths,nz,*rowstart,itmp;
1301be0148SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0;
14dd6ea824SBarry Smith   MatScalar      *av,*bv;
154e5e7fe4SHong Zhang 
164e5e7fe4SHong Zhang   PetscFunctionBegin;
174e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
18dcca6d9dSJed Brown   ierr = PetscMalloc2(m,&rowlengths,m+1,&rowstart);CHKERRQ(ierr);
19a7a3a9ebSHong Zhang 
20a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
21a7a3a9ebSHong Zhang   k  = 0;
22a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
234e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
2401be0148SBarry Smith     if (nz) {
2501be0148SBarry Smith       rowlengths[k] += nz;   /* no. of upper triangular blocks */
2601be0148SBarry Smith       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
2701be0148SBarry Smith       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
28a7a3a9ebSHong Zhang         rowlengths[(*aj)*bs]++; aj++;
29a7a3a9ebSHong Zhang       }
3001be0148SBarry Smith     }
31a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
32a7a3a9ebSHong Zhang     for (j=1; j<bs; j++) {
33a7a3a9ebSHong Zhang       rowlengths[k+j] = rowlengths[k];
34a7a3a9ebSHong Zhang     }
35a7a3a9ebSHong Zhang     k += bs;
364e5e7fe4SHong Zhang   }
374e5e7fe4SHong Zhang 
38*bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
39ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
40f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
416d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
42f204ca49SKris Buschelman     ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
43*bd019fc1SStefano Zampini     ierr = MatSetBlockSize(B,A->rmap->bs);CHKERRQ(ierr);
44*bd019fc1SStefano Zampini   } else B = *newmat;
454e5e7fe4SHong Zhang 
464e5e7fe4SHong Zhang   b  = (Mat_SeqAIJ*)(B->data);
474e5e7fe4SHong Zhang   bi = b->i;
484e5e7fe4SHong Zhang   bj = b->j;
494e5e7fe4SHong Zhang   bv = b->a;
504e5e7fe4SHong Zhang 
514e5e7fe4SHong Zhang   /* set b->i */
52a7a3a9ebSHong Zhang   bi[0] = 0; rowstart[0] = 0;
53a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
54a7a3a9ebSHong Zhang     for (j=0; j<bs; j++) {
55a7a3a9ebSHong Zhang       b->ilen[i*bs+j]    = rowlengths[i*bs];
56a7a3a9ebSHong Zhang       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
574e5e7fe4SHong Zhang     }
58a7a3a9ebSHong Zhang     bi[i+1] = bi[i] + rowlengths[i*bs]/bs;
59a7a3a9ebSHong Zhang   }
6001be0148SBarry Smith   if (bi[mbs] != 2*a->nz - diagcnt) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-diagcnt: %D\n",bi[mbs],2*a->nz - diagcnt);
614e5e7fe4SHong Zhang 
624e5e7fe4SHong Zhang   /* set b->j and b->a */
634e5e7fe4SHong Zhang   aj = a->j; av = a->a;
64a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
6501be0148SBarry Smith     nz = ai[i+1] - ai[i];
66a7a3a9ebSHong Zhang     /* diagonal block */
6701be0148SBarry Smith     if (nz && *aj == i) {
6801be0148SBarry Smith       nz--;
69a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row i*bs+j */
70a7a3a9ebSHong Zhang         itmp = i*bs+j;
71a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col i*bs+k */
72a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
73a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
74a7a3a9ebSHong Zhang           rowstart[itmp]++;
75a7a3a9ebSHong Zhang         }
76a7a3a9ebSHong Zhang       }
77a7a3a9ebSHong Zhang       aj++; av += bs2;
7801be0148SBarry Smith     }
79a7a3a9ebSHong Zhang 
804e5e7fe4SHong Zhang     while (nz--) {
81a7a3a9ebSHong Zhang       /* lower triangular blocks */
82a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row (*aj)*bs+j */
83a7a3a9ebSHong Zhang         itmp = (*aj)*bs+j;
84a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col i*bs+k */
85a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = i*bs+k;
863bededecSBarry Smith           *(bv + rowstart[itmp]) = *(av+j*bs+k);
87a7a3a9ebSHong Zhang           rowstart[itmp]++;
88a7a3a9ebSHong Zhang         }
89a7a3a9ebSHong Zhang       }
90a7a3a9ebSHong Zhang       /* upper triangular blocks */
91a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row i*bs+j */
92a7a3a9ebSHong Zhang         itmp = i*bs+j;
93a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col (*aj)*bs+k */
94a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
95a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
96a7a3a9ebSHong Zhang           rowstart[itmp]++;
97a7a3a9ebSHong Zhang         }
98a7a3a9ebSHong Zhang       }
99a7a3a9ebSHong Zhang       aj++; av += bs2;
1004e5e7fe4SHong Zhang     }
1014e5e7fe4SHong Zhang   }
10274ed9c26SBarry Smith   ierr = PetscFree2(rowlengths,rowstart);CHKERRQ(ierr);
1034e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1044e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1054e5e7fe4SHong Zhang 
106511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
10728be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
108c3d102feSKris Buschelman   } else {
1094e5e7fe4SHong Zhang     *newmat = B;
110c3d102feSKris Buschelman   }
1114e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1124e5e7fe4SHong Zhang }
113be1d678aSKris Buschelman 
114cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1150adebc6cSBarry Smith {
116676c34cdSKris Buschelman   Mat            B;
11759557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
118861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
119dfbe8321SBarry Smith   PetscErrorCode ierr;
120d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
121dd6ea824SBarry Smith   MatScalar      *av,*bv;
12259557b74SHong Zhang 
12359557b74SHong Zhang   PetscFunctionBegin;
124ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
125e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
12659557b74SHong Zhang 
127785e854fSJed Brown   ierr = PetscMalloc1(m,&rowlengths);CHKERRQ(ierr);
12859557b74SHong Zhang   for (i=0; i<m; i++) {
12959557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13059557b74SHong Zhang   }
131*bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
132ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
133f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1346d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
135367daffbSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr);
136*bd019fc1SStefano Zampini   } else B = *newmat;
13759557b74SHong Zhang 
138676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
139861ba921SHong Zhang   bi = b->i;
140861ba921SHong Zhang   bj = b->j;
141861ba921SHong Zhang   bv = b->a;
142861ba921SHong Zhang 
143861ba921SHong Zhang   bi[0] = 0;
14459557b74SHong Zhang   for (i=0; i<m; i++) {
14559557b74SHong Zhang     aj = a->j + a->diag[i];
14659557b74SHong Zhang     av = a->a + a->diag[i];
147861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++) {
148861ba921SHong Zhang       *bj = *aj; bj++; aj++;
149861ba921SHong Zhang       *bv = *av; bv++; av++;
150861ba921SHong Zhang     }
151861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
152861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
15359557b74SHong Zhang   }
15459557b74SHong Zhang 
15559557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
156676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158676c34cdSKris Buschelman 
159511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16028be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
161c3d102feSKris Buschelman   } else {
162676c34cdSKris Buschelman     *newmat = B;
163c3d102feSKris Buschelman   }
16459557b74SHong Zhang   PetscFunctionReturn(0);
16559557b74SHong Zhang }
16659557b74SHong Zhang 
167cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
168a0e1a404SHong Zhang {
169a0e1a404SHong Zhang   Mat            B;
170a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
171a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
172dfbe8321SBarry Smith   PetscErrorCode ierr;
173d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
17474ee4d9fSHong Zhang   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
175dd6ea824SBarry Smith   MatScalar      *av,*bv;
176a0e1a404SHong Zhang 
177a0e1a404SHong Zhang   PetscFunctionBegin;
178a0e1a404SHong Zhang   /* compute browlengths of newmat */
179dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&browlengths,mbs,&browstart);CHKERRQ(ierr);
180a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
181a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
182a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
183a0e1a404SHong Zhang     aj++; /* skip diagonal */
184a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
185a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
186a0e1a404SHong Zhang     }
187a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
188a0e1a404SHong Zhang   }
189a0e1a404SHong Zhang 
190*bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
191ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
192f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1936d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
194f204ca49SKris Buschelman     ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
195*bd019fc1SStefano Zampini   } else B = *newmat;
196a0e1a404SHong Zhang 
197a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
198a0e1a404SHong Zhang   bi = b->i;
199a0e1a404SHong Zhang   bj = b->j;
200a0e1a404SHong Zhang   bv = b->a;
201a0e1a404SHong Zhang 
202a0e1a404SHong Zhang   /* set b->i */
203a0e1a404SHong Zhang   bi[0] = 0;
204a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
205a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
206a0e1a404SHong Zhang     bi[i+1]      = bi[i] + browlengths[i];
207a0e1a404SHong Zhang     browstart[i] = bi[i];
208a0e1a404SHong Zhang   }
209e32f2f54SBarry Smith   if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs);
210a0e1a404SHong Zhang 
211a0e1a404SHong Zhang   /* set b->j and b->a */
212a0e1a404SHong Zhang   aj = a->j; av = a->a;
213a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
214a0e1a404SHong Zhang     /* diagonal block */
215a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
21626fbe8dcSKarl Rupp 
217a0e1a404SHong Zhang     itmp = bs2*browstart[i];
218a0e1a404SHong Zhang     for (k=0; k<bs2; k++) {
219a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
220a0e1a404SHong Zhang     }
221a0e1a404SHong Zhang     browstart[i]++;
222a0e1a404SHong Zhang 
223a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
224a0e1a404SHong Zhang     while (nz--) {
22574ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
22674ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
22726fbe8dcSKarl Rupp 
22874ee4d9fSHong Zhang       itmp = bs2*browstart[*aj];  /* row index */
22974ee4d9fSHong Zhang       for (col=0; col<bs; col++) {
23074ee4d9fSHong Zhang         k = col;
23174ee4d9fSHong Zhang         for (row=0; row<bs; row++) {
23274ee4d9fSHong Zhang           bv[itmp + col*bs+row] = av[k]; k+=bs;
23374ee4d9fSHong Zhang         }
234a0e1a404SHong Zhang       }
235a0e1a404SHong Zhang       browstart[*aj]++;
236a0e1a404SHong Zhang 
237a0e1a404SHong Zhang       /* upper triangular blocks */
238a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
23926fbe8dcSKarl Rupp 
240a0e1a404SHong Zhang       itmp = bs2*browstart[i];
241a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
24274ee4d9fSHong Zhang         bv[itmp + k] = av[k];
243a0e1a404SHong Zhang       }
24474ee4d9fSHong Zhang       av += bs2;
245a0e1a404SHong Zhang       browstart[i]++;
246a0e1a404SHong Zhang     }
247a0e1a404SHong Zhang   }
24874ed9c26SBarry Smith   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
249a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
250a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251a0e1a404SHong Zhang 
252511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
25328be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
254c3d102feSKris Buschelman   } else {
255a0e1a404SHong Zhang     *newmat = B;
256c3d102feSKris Buschelman   }
257a0e1a404SHong Zhang   PetscFunctionReturn(0);
258a0e1a404SHong Zhang }
259be1d678aSKris Buschelman 
260cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
261a0e1a404SHong Zhang {
262a0e1a404SHong Zhang   Mat            B;
263a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
264a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
265dfbe8321SBarry Smith   PetscErrorCode ierr;
266d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
267d0f46423SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
268dd6ea824SBarry Smith   MatScalar      *av,*bv;
269ace3abfcSBarry Smith   PetscBool      flg;
270a0e1a404SHong Zhang 
271a0e1a404SHong Zhang   PetscFunctionBegin;
272ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
273e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2742af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
275e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
276a0e1a404SHong Zhang 
277785e854fSJed Brown   ierr = PetscMalloc1(mbs,&browlengths);CHKERRQ(ierr);
278a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
279a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
280a0e1a404SHong Zhang   }
281a0e1a404SHong Zhang 
282*bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
283ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
284f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2856d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
286367daffbSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
287*bd019fc1SStefano Zampini   } else B = *newmat;
288a0e1a404SHong Zhang 
289a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
290a0e1a404SHong Zhang   bi = b->i;
291a0e1a404SHong Zhang   bj = b->j;
292a0e1a404SHong Zhang   bv = b->a;
293a0e1a404SHong Zhang 
294a0e1a404SHong Zhang   bi[0] = 0;
295a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
296a0e1a404SHong Zhang     aj = a->j + a->diag[i];
297a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
298a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++) {
299a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
300a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
301a0e1a404SHong Zhang         *bv = *av; bv++; av++;
302a0e1a404SHong Zhang       }
303a0e1a404SHong Zhang     }
304a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
305a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
306a0e1a404SHong Zhang   }
307a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
308a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
309a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
310a0e1a404SHong Zhang 
311511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
31228be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
313c3d102feSKris Buschelman   } else {
314a0e1a404SHong Zhang     *newmat = B;
315c3d102feSKris Buschelman   }
316a0e1a404SHong Zhang   PetscFunctionReturn(0);
317a0e1a404SHong Zhang }
318