xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
259557b74SHong Zhang 
3*7c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
4*7c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
5*7c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
659557b74SHong Zhang 
759557b74SHong Zhang EXTERN_C_BEGIN
859557b74SHong Zhang #undef __FUNCT__
94e5e7fe4SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqAIJ"
10f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
114e5e7fe4SHong Zhang {
124e5e7fe4SHong Zhang   Mat            B;
134e5e7fe4SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
144e5e7fe4SHong Zhang   Mat_SeqAIJ     *b;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
16d0f46423SBarry 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;
17d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs;
18dd6ea824SBarry Smith   MatScalar      *av,*bv;
194e5e7fe4SHong Zhang 
204e5e7fe4SHong Zhang   PetscFunctionBegin;
214e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
2213f74950SBarry Smith   ierr = PetscMalloc((2*m+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
23a7a3a9ebSHong Zhang   rowstart = rowlengths + m;
24a7a3a9ebSHong Zhang 
25a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
264e5e7fe4SHong Zhang   aj = a->j;
27a7a3a9ebSHong Zhang   k = 0;
28a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
294e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
30a7a3a9ebSHong Zhang     aj++; /* skip diagonal */
31a7a3a9ebSHong Zhang     for (j=1; j<nz; j++) { /* no. of lower triangular blocks */
32a7a3a9ebSHong Zhang       rowlengths[(*aj)*bs]++; aj++;
33a7a3a9ebSHong Zhang     }
34a7a3a9ebSHong Zhang     rowlengths[k] += nz;   /* no. of upper triangular blocks */
35a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
36a7a3a9ebSHong Zhang     for (j=1; j<bs; j++) {
37a7a3a9ebSHong Zhang       rowlengths[k+j] = rowlengths[k];
38a7a3a9ebSHong Zhang     }
39a7a3a9ebSHong Zhang     k += bs;
40a7a3a9ebSHong Zhang     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
414e5e7fe4SHong Zhang   }
424e5e7fe4SHong Zhang 
437adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
44f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
45f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
46f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
474e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
48d0f46423SBarry Smith   B->rmap->bs = A->rmap->bs;
494e5e7fe4SHong Zhang 
504e5e7fe4SHong Zhang   b  = (Mat_SeqAIJ*)(B->data);
514e5e7fe4SHong Zhang   bi = b->i;
524e5e7fe4SHong Zhang   bj = b->j;
534e5e7fe4SHong Zhang   bv = b->a;
544e5e7fe4SHong Zhang 
554e5e7fe4SHong Zhang   /* set b->i */
56a7a3a9ebSHong Zhang   bi[0] = 0; rowstart[0] = 0;
57a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++){
58a7a3a9ebSHong Zhang     for (j=0; j<bs; j++){
59a7a3a9ebSHong Zhang       b->ilen[i*bs+j]  = rowlengths[i*bs];
60a7a3a9ebSHong Zhang       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
614e5e7fe4SHong Zhang     }
62a7a3a9ebSHong Zhang     bi[i+1]     = bi[i] + rowlengths[i*bs]/bs;
63a7a3a9ebSHong Zhang   }
6477431f27SBarry Smith   if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-mbs: %D\n",bi[mbs],2*a->nz - mbs);
654e5e7fe4SHong Zhang 
664e5e7fe4SHong Zhang   /* set b->j and b->a */
674e5e7fe4SHong Zhang   aj = a->j; av = a->a;
68a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
69a7a3a9ebSHong Zhang     /* diagonal block */
70a7a3a9ebSHong Zhang     for (j=0; j<bs; j++){   /* row i*bs+j */
71a7a3a9ebSHong Zhang       itmp = i*bs+j;
72a7a3a9ebSHong Zhang       for (k=0; k<bs; k++){ /* col i*bs+k */
73a7a3a9ebSHong Zhang         *(bj + rowstart[itmp]) = (*aj)*bs+k;
74a7a3a9ebSHong Zhang         *(bv + rowstart[itmp]) = *(av+k*bs+j);
75a7a3a9ebSHong Zhang         rowstart[itmp]++;
76a7a3a9ebSHong Zhang       }
77a7a3a9ebSHong Zhang     }
78a7a3a9ebSHong Zhang     aj++; av += bs2;
79a7a3a9ebSHong Zhang 
804e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i] -1;
814e5e7fe4SHong Zhang     while (nz--){
82a7a3a9ebSHong Zhang       /* lower triangular blocks */
83a7a3a9ebSHong Zhang       for (j=0; j<bs; j++){   /* row (*aj)*bs+j */
84a7a3a9ebSHong Zhang         itmp = (*aj)*bs+j;
85a7a3a9ebSHong Zhang         for (k=0; k<bs; k++){ /* col i*bs+k */
86a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = i*bs+k;
87a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
88a7a3a9ebSHong Zhang           rowstart[itmp]++;
89a7a3a9ebSHong Zhang         }
90a7a3a9ebSHong Zhang       }
91a7a3a9ebSHong Zhang       /* upper triangular blocks */
92a7a3a9ebSHong Zhang       for (j=0; j<bs; j++){   /* row i*bs+j */
93a7a3a9ebSHong Zhang         itmp = i*bs+j;
94a7a3a9ebSHong Zhang         for (k=0; k<bs; k++){ /* col (*aj)*bs+k */
95a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
96a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
97a7a3a9ebSHong Zhang           rowstart[itmp]++;
98a7a3a9ebSHong Zhang         }
99a7a3a9ebSHong Zhang       }
100a7a3a9ebSHong Zhang       aj++; av += bs2;
1014e5e7fe4SHong Zhang     }
1024e5e7fe4SHong Zhang   }
1034e5e7fe4SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1044e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1054e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1064e5e7fe4SHong Zhang 
107ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1088ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
109c3d102feSKris Buschelman   } else {
1104e5e7fe4SHong Zhang     *newmat = B;
111c3d102feSKris Buschelman   }
1124e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1134e5e7fe4SHong Zhang }
114be1d678aSKris Buschelman EXTERN_C_END
115be1d678aSKris Buschelman 
116be1d678aSKris Buschelman EXTERN_C_BEGIN
1174e5e7fe4SHong Zhang #undef __FUNCT__
11859557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
1198cf70c4bSSatish Balay PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) {
120676c34cdSKris Buschelman   Mat            B;
12159557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
122861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
123dfbe8321SBarry Smith   PetscErrorCode ierr;
124d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
125dd6ea824SBarry Smith   MatScalar      *av,*bv;
12659557b74SHong Zhang 
12759557b74SHong Zhang   PetscFunctionBegin;
1282d9a3abdSHong Zhang   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
12959557b74SHong Zhang 
13013f74950SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13159557b74SHong Zhang   for (i=0; i<m; i++) {
13259557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13359557b74SHong Zhang   }
1347adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
135f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
136f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
137ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
13859557b74SHong Zhang 
1394e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
14059557b74SHong Zhang 
141676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
142861ba921SHong Zhang   bi = b->i;
143861ba921SHong Zhang   bj = b->j;
144861ba921SHong Zhang   bv = b->a;
145861ba921SHong Zhang 
146861ba921SHong Zhang   bi[0] = 0;
14759557b74SHong Zhang   for (i=0; i<m; i++) {
14859557b74SHong Zhang     aj = a->j + a->diag[i];
14959557b74SHong Zhang     av = a->a + a->diag[i];
150861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++){
151861ba921SHong Zhang       *bj = *aj; bj++; aj++;
152861ba921SHong Zhang       *bv = *av; bv++; av++;
153861ba921SHong Zhang     }
154861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
155861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
15659557b74SHong Zhang   }
15759557b74SHong Zhang 
15859557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
159676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161676c34cdSKris Buschelman 
162ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1638ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
164c3d102feSKris Buschelman   } else {
165676c34cdSKris Buschelman     *newmat = B;
166c3d102feSKris Buschelman   }
16759557b74SHong Zhang   PetscFunctionReturn(0);
16859557b74SHong Zhang }
16959557b74SHong Zhang EXTERN_C_END
17059557b74SHong Zhang 
171a0e1a404SHong Zhang EXTERN_C_BEGIN
172a0e1a404SHong Zhang #undef __FUNCT__
173a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ"
174f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
175a0e1a404SHong Zhang {
176a0e1a404SHong Zhang   Mat            B;
177a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
178a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
179dfbe8321SBarry Smith   PetscErrorCode ierr;
180d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
181d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs;
182dd6ea824SBarry Smith   MatScalar      *av,*bv;
183a0e1a404SHong Zhang 
184a0e1a404SHong Zhang   PetscFunctionBegin;
185a0e1a404SHong Zhang   /* compute browlengths of newmat */
18613f74950SBarry Smith   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
187a0e1a404SHong Zhang   browstart = browlengths + mbs;
188a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
189a0e1a404SHong Zhang   aj = a->j;
190a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
191a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
192a0e1a404SHong Zhang     aj++; /* skip diagonal */
193a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
194a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
195a0e1a404SHong Zhang     }
196a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
197a0e1a404SHong Zhang   }
198a0e1a404SHong Zhang 
1997adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
200f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
201f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
202f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
2034e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
204a0e1a404SHong Zhang 
205a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
206a0e1a404SHong Zhang   bi = b->i;
207a0e1a404SHong Zhang   bj = b->j;
208a0e1a404SHong Zhang   bv = b->a;
209a0e1a404SHong Zhang 
210a0e1a404SHong Zhang   /* set b->i */
211a0e1a404SHong Zhang   bi[0] = 0;
212a0e1a404SHong Zhang   for (i=0; i<mbs; i++){
213a0e1a404SHong Zhang     b->ilen[i]    = browlengths[i];
214a0e1a404SHong Zhang     bi[i+1]       = bi[i] + browlengths[i];
215a0e1a404SHong Zhang     browstart[i]  = bi[i];
216a0e1a404SHong Zhang   }
21777431f27SBarry Smith   if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs);
218a0e1a404SHong Zhang 
219a0e1a404SHong Zhang   /* set b->j and b->a */
220a0e1a404SHong Zhang   aj = a->j; av = a->a;
221a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
222a0e1a404SHong Zhang     /* diagonal block */
223a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
224a0e1a404SHong Zhang     itmp = bs2*browstart[i];
225a0e1a404SHong Zhang     for (k=0; k<bs2; k++){
226a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
227a0e1a404SHong Zhang     }
228a0e1a404SHong Zhang     browstart[i]++;
229a0e1a404SHong Zhang 
230a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
231a0e1a404SHong Zhang     while (nz--){
232a0e1a404SHong Zhang       /* lower triangular blocks */
233a0e1a404SHong Zhang       *(bj + browstart[*aj]) = i;
234a0e1a404SHong Zhang       itmp = bs2*browstart[*aj];
235a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
236a0e1a404SHong Zhang         *(bv + itmp + k) = *(av + k);
237a0e1a404SHong Zhang       }
238a0e1a404SHong Zhang       browstart[*aj]++;
239a0e1a404SHong Zhang 
240a0e1a404SHong Zhang       /* upper triangular blocks */
241a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
242a0e1a404SHong Zhang       itmp = bs2*browstart[i];
243a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
244a0e1a404SHong Zhang         *(bv + itmp + k) = *av; av++;
245a0e1a404SHong Zhang       }
246a0e1a404SHong Zhang       browstart[i]++;
247a0e1a404SHong Zhang     }
248a0e1a404SHong Zhang   }
249a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
250a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252a0e1a404SHong Zhang 
253ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
2548ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
255c3d102feSKris Buschelman   } else {
256a0e1a404SHong Zhang     *newmat = B;
257c3d102feSKris Buschelman   }
258a0e1a404SHong Zhang   PetscFunctionReturn(0);
259a0e1a404SHong Zhang }
260be1d678aSKris Buschelman EXTERN_C_END
261be1d678aSKris Buschelman 
262be1d678aSKris Buschelman EXTERN_C_BEGIN
263a0e1a404SHong Zhang #undef __FUNCT__
264a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
265f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
266a0e1a404SHong Zhang {
267a0e1a404SHong Zhang   Mat            B;
268a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
269a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
270dfbe8321SBarry Smith   PetscErrorCode ierr;
271d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
272d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
273dd6ea824SBarry Smith   MatScalar      *av,*bv;
2742af78befSBarry Smith   PetscTruth     flg;
275a0e1a404SHong Zhang 
276a0e1a404SHong Zhang   PetscFunctionBegin;
277a0e1a404SHong Zhang   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
2782af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
2792af78befSBarry Smith   if (flg) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
280a0e1a404SHong Zhang 
28113f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
282a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
283a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
284a0e1a404SHong Zhang   }
285a0e1a404SHong Zhang 
2867adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
287f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
288f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
289ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr);
2904e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
291a0e1a404SHong Zhang 
292a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
293a0e1a404SHong Zhang   bi = b->i;
294a0e1a404SHong Zhang   bj = b->j;
295a0e1a404SHong Zhang   bv = b->a;
296a0e1a404SHong Zhang 
297a0e1a404SHong Zhang   bi[0] = 0;
298a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
299a0e1a404SHong Zhang     aj = a->j + a->diag[i];
300a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
301a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++){
302a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
303a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
304a0e1a404SHong Zhang         *bv = *av; bv++; av++;
305a0e1a404SHong Zhang       }
306a0e1a404SHong Zhang     }
307a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
308a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
309a0e1a404SHong Zhang   }
310a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
311a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313a0e1a404SHong Zhang 
314ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
3158ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
316c3d102feSKris Buschelman   } else {
317a0e1a404SHong Zhang     *newmat = B;
318c3d102feSKris Buschelman   }
319a0e1a404SHong Zhang 
320a0e1a404SHong Zhang   PetscFunctionReturn(0);
321a0e1a404SHong Zhang }
322a0e1a404SHong Zhang EXTERN_C_END
323