xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
259557b74SHong Zhang 
359557b74SHong Zhang #include "src/mat/impls/aij/seq/aij.h"
4595208bcSHong Zhang #include "src/mat/impls/baij/seq/baij.h"
5861ba921SHong Zhang #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"
10*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat A,const 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;
1613f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,j,k,*bi,*bj,*rowlengths,nz,*rowstart,itmp;
1713f74950SBarry Smith   PetscInt       bs=A->bs,bs2=bs*bs,mbs=A->m/bs;
184e5e7fe4SHong Zhang   PetscScalar    *av,*bv;
194e5e7fe4SHong Zhang 
204e5e7fe4SHong Zhang   PetscFunctionBegin;
21a7a3a9ebSHong Zhang 
224e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
2313f74950SBarry Smith   ierr = PetscMalloc((2*m+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
24a7a3a9ebSHong Zhang   rowstart = rowlengths + m;
25a7a3a9ebSHong Zhang 
26a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
274e5e7fe4SHong Zhang   aj = a->j;
28a7a3a9ebSHong Zhang   k = 0;
29a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
304e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
31a7a3a9ebSHong Zhang     aj++; /* skip diagonal */
32a7a3a9ebSHong Zhang     for (j=1; j<nz; j++) { /* no. of lower triangular blocks */
33a7a3a9ebSHong Zhang       rowlengths[(*aj)*bs]++; aj++;
34a7a3a9ebSHong Zhang     }
35a7a3a9ebSHong Zhang     rowlengths[k] += nz;   /* no. of upper triangular blocks */
36a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
37a7a3a9ebSHong Zhang     for (j=1; j<bs; j++) {
38a7a3a9ebSHong Zhang       rowlengths[k+j] = rowlengths[k];
39a7a3a9ebSHong Zhang     }
40a7a3a9ebSHong Zhang     k += bs;
41a7a3a9ebSHong Zhang     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
424e5e7fe4SHong Zhang   }
434e5e7fe4SHong Zhang 
44f204ca49SKris Buschelman   ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr);
45f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
46f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
47a7a3a9ebSHong Zhang   ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
484e5e7fe4SHong Zhang   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
494e5e7fe4SHong Zhang   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
50521d7252SBarry Smith   B->bs = A->bs;
514e5e7fe4SHong Zhang 
524e5e7fe4SHong Zhang   b  = (Mat_SeqAIJ*)(B->data);
534e5e7fe4SHong Zhang   bi = b->i;
544e5e7fe4SHong Zhang   bj = b->j;
554e5e7fe4SHong Zhang   bv = b->a;
564e5e7fe4SHong Zhang 
574e5e7fe4SHong Zhang   /* set b->i */
58a7a3a9ebSHong Zhang   bi[0] = 0; rowstart[0] = 0;
59a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++){
60a7a3a9ebSHong Zhang     for (j=0; j<bs; j++){
61a7a3a9ebSHong Zhang       b->ilen[i*bs+j]  = rowlengths[i*bs];
62a7a3a9ebSHong Zhang       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
634e5e7fe4SHong Zhang     }
64a7a3a9ebSHong Zhang     bi[i+1]     = bi[i] + rowlengths[i*bs]/bs;
65a7a3a9ebSHong Zhang   }
6677431f27SBarry 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);
674e5e7fe4SHong Zhang 
684e5e7fe4SHong Zhang   /* set b->j and b->a */
694e5e7fe4SHong Zhang   aj = a->j; av = a->a;
70a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
71a7a3a9ebSHong Zhang     /* diagonal block */
72a7a3a9ebSHong Zhang     for (j=0; j<bs; j++){   /* row i*bs+j */
73a7a3a9ebSHong Zhang       itmp = i*bs+j;
74a7a3a9ebSHong Zhang       for (k=0; k<bs; k++){ /* col i*bs+k */
75a7a3a9ebSHong Zhang         *(bj + rowstart[itmp]) = (*aj)*bs+k;
76a7a3a9ebSHong Zhang         *(bv + rowstart[itmp]) = *(av+k*bs+j);
77a7a3a9ebSHong Zhang         rowstart[itmp]++;
78a7a3a9ebSHong Zhang       }
79a7a3a9ebSHong Zhang     }
80a7a3a9ebSHong Zhang     aj++; av += bs2;
81a7a3a9ebSHong Zhang 
824e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i] -1;
834e5e7fe4SHong Zhang     while (nz--){
84a7a3a9ebSHong Zhang       /* lower triangular blocks */
85a7a3a9ebSHong Zhang       for (j=0; j<bs; j++){   /* row (*aj)*bs+j */
86a7a3a9ebSHong Zhang         itmp = (*aj)*bs+j;
87a7a3a9ebSHong Zhang         for (k=0; k<bs; k++){ /* col i*bs+k */
88a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = i*bs+k;
89a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
90a7a3a9ebSHong Zhang           rowstart[itmp]++;
91a7a3a9ebSHong Zhang         }
92a7a3a9ebSHong Zhang       }
93a7a3a9ebSHong Zhang       /* upper triangular blocks */
94a7a3a9ebSHong Zhang       for (j=0; j<bs; j++){   /* row i*bs+j */
95a7a3a9ebSHong Zhang         itmp = i*bs+j;
96a7a3a9ebSHong Zhang         for (k=0; k<bs; k++){ /* col (*aj)*bs+k */
97a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
98a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
99a7a3a9ebSHong Zhang           rowstart[itmp]++;
100a7a3a9ebSHong Zhang         }
101a7a3a9ebSHong Zhang       }
102a7a3a9ebSHong Zhang       aj++; av += bs2;
1034e5e7fe4SHong Zhang     }
1044e5e7fe4SHong Zhang   }
1054e5e7fe4SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1064e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1074e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1084e5e7fe4SHong Zhang 
109ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1108ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
111c3d102feSKris Buschelman   } else {
1124e5e7fe4SHong Zhang     *newmat = B;
113c3d102feSKris Buschelman   }
1144e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1154e5e7fe4SHong Zhang }
116*be1d678aSKris Buschelman EXTERN_C_END
117*be1d678aSKris Buschelman 
118*be1d678aSKris Buschelman EXTERN_C_BEGIN
1194e5e7fe4SHong Zhang #undef __FUNCT__
12059557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
121*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) {
122676c34cdSKris Buschelman   Mat            B;
12359557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
124861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
125dfbe8321SBarry Smith   PetscErrorCode ierr;
12679416396SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->M,n=A->N,i,j,*bi,*bj,*rowlengths;
127861ba921SHong Zhang   PetscScalar    *av,*bv;
12859557b74SHong Zhang 
12959557b74SHong Zhang   PetscFunctionBegin;
1302d9a3abdSHong Zhang   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
13118f87311SHong Zhang   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
13259557b74SHong Zhang 
13313f74950SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13459557b74SHong Zhang   for (i=0; i<m; i++) {
13559557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13659557b74SHong Zhang   }
137f204ca49SKris Buschelman   ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr);
138f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
139f204ca49SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr);
14059557b74SHong Zhang 
141676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
142676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
143676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
14459557b74SHong Zhang 
145676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
146861ba921SHong Zhang   bi = b->i;
147861ba921SHong Zhang   bj = b->j;
148861ba921SHong Zhang   bv = b->a;
149861ba921SHong Zhang 
150861ba921SHong Zhang   bi[0] = 0;
15159557b74SHong Zhang   for (i=0; i<m; i++) {
15259557b74SHong Zhang     aj = a->j + a->diag[i];
15359557b74SHong Zhang     av = a->a + a->diag[i];
154861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++){
155861ba921SHong Zhang       *bj = *aj; bj++; aj++;
156861ba921SHong Zhang       *bv = *av; bv++; av++;
157861ba921SHong Zhang     }
158861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
159861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
16059557b74SHong Zhang   }
16159557b74SHong Zhang 
16259557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
163676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165676c34cdSKris Buschelman 
166ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1678ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
168c3d102feSKris Buschelman   } else {
169676c34cdSKris Buschelman     *newmat = B;
170c3d102feSKris Buschelman   }
17159557b74SHong Zhang   PetscFunctionReturn(0);
17259557b74SHong Zhang }
17359557b74SHong Zhang EXTERN_C_END
17459557b74SHong Zhang 
175a0e1a404SHong Zhang EXTERN_C_BEGIN
176a0e1a404SHong Zhang #undef __FUNCT__
177a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ"
178*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
179a0e1a404SHong Zhang {
180a0e1a404SHong Zhang   Mat            B;
181a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
182a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
183dfbe8321SBarry Smith   PetscErrorCode ierr;
18413f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
18513f74950SBarry Smith   PetscInt       bs=A->bs,bs2=bs*bs,mbs=m/bs;
186a0e1a404SHong Zhang   PetscScalar    *av,*bv;
187a0e1a404SHong Zhang 
188a0e1a404SHong Zhang   PetscFunctionBegin;
189a0e1a404SHong Zhang   /* compute browlengths of newmat */
19013f74950SBarry Smith   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
191a0e1a404SHong Zhang   browstart = browlengths + mbs;
192a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
193a0e1a404SHong Zhang   aj = a->j;
194a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
195a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
196a0e1a404SHong Zhang     aj++; /* skip diagonal */
197a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
198a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
199a0e1a404SHong Zhang     }
200a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
201a0e1a404SHong Zhang   }
202a0e1a404SHong Zhang 
203f204ca49SKris Buschelman   ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr);
204f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
205f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
206a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
207a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
208a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
209a0e1a404SHong Zhang 
210a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
211a0e1a404SHong Zhang   bi = b->i;
212a0e1a404SHong Zhang   bj = b->j;
213a0e1a404SHong Zhang   bv = b->a;
214a0e1a404SHong Zhang 
215a0e1a404SHong Zhang   /* set b->i */
216a0e1a404SHong Zhang   bi[0] = 0;
217a0e1a404SHong Zhang   for (i=0; i<mbs; i++){
218a0e1a404SHong Zhang     b->ilen[i]    = browlengths[i];
219a0e1a404SHong Zhang     bi[i+1]       = bi[i] + browlengths[i];
220a0e1a404SHong Zhang     browstart[i]  = bi[i];
221a0e1a404SHong Zhang   }
22277431f27SBarry 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);
223a0e1a404SHong Zhang 
224a0e1a404SHong Zhang   /* set b->j and b->a */
225a0e1a404SHong Zhang   aj = a->j; av = a->a;
226a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
227a0e1a404SHong Zhang     /* diagonal block */
228a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
229a0e1a404SHong Zhang     itmp = bs2*browstart[i];
230a0e1a404SHong Zhang     for (k=0; k<bs2; k++){
231a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
232a0e1a404SHong Zhang     }
233a0e1a404SHong Zhang     browstart[i]++;
234a0e1a404SHong Zhang 
235a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
236a0e1a404SHong Zhang     while (nz--){
237a0e1a404SHong Zhang       /* lower triangular blocks */
238a0e1a404SHong Zhang       *(bj + browstart[*aj]) = i;
239a0e1a404SHong Zhang       itmp = bs2*browstart[*aj];
240a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
241a0e1a404SHong Zhang         *(bv + itmp + k) = *(av + k);
242a0e1a404SHong Zhang       }
243a0e1a404SHong Zhang       browstart[*aj]++;
244a0e1a404SHong Zhang 
245a0e1a404SHong Zhang       /* upper triangular blocks */
246a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
247a0e1a404SHong Zhang       itmp = bs2*browstart[i];
248a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
249a0e1a404SHong Zhang         *(bv + itmp + k) = *av; av++;
250a0e1a404SHong Zhang       }
251a0e1a404SHong Zhang       browstart[i]++;
252a0e1a404SHong Zhang     }
253a0e1a404SHong Zhang   }
254a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
255a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257a0e1a404SHong Zhang 
258ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
2598ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
260c3d102feSKris Buschelman   } else {
261a0e1a404SHong Zhang     *newmat = B;
262c3d102feSKris Buschelman   }
263a0e1a404SHong Zhang   PetscFunctionReturn(0);
264a0e1a404SHong Zhang }
265*be1d678aSKris Buschelman EXTERN_C_END
266*be1d678aSKris Buschelman 
267*be1d678aSKris Buschelman EXTERN_C_BEGIN
268a0e1a404SHong Zhang #undef __FUNCT__
269a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
270*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqBAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
271a0e1a404SHong Zhang {
272a0e1a404SHong Zhang   Mat            B;
273a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
274a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
275dfbe8321SBarry Smith   PetscErrorCode ierr;
27613f74950SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->m,n=A->n,i,j,k,*bi,*bj,*browlengths;
27713f74950SBarry Smith   PetscInt       bs=A->bs,bs2=bs*bs,mbs=m/bs;
278a0e1a404SHong Zhang   PetscScalar    *av,*bv;
279a0e1a404SHong Zhang 
280a0e1a404SHong Zhang   PetscFunctionBegin;
281a0e1a404SHong Zhang   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
282595208bcSHong Zhang   ierr = MatMissingDiagonal_SeqBAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
283a0e1a404SHong Zhang 
28413f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
285a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
286a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
287a0e1a404SHong Zhang   }
288a0e1a404SHong Zhang 
289f204ca49SKris Buschelman   ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr);
290f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
291f204ca49SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
292a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
293a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
294a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
295a0e1a404SHong Zhang 
296a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
297a0e1a404SHong Zhang   bi = b->i;
298a0e1a404SHong Zhang   bj = b->j;
299a0e1a404SHong Zhang   bv = b->a;
300a0e1a404SHong Zhang 
301a0e1a404SHong Zhang   bi[0] = 0;
302a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
303a0e1a404SHong Zhang     aj = a->j + a->diag[i];
304a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
305a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++){
306a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
307a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
308a0e1a404SHong Zhang         *bv = *av; bv++; av++;
309a0e1a404SHong Zhang       }
310a0e1a404SHong Zhang     }
311a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
312a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
313a0e1a404SHong Zhang   }
314a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
315a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
316a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
317a0e1a404SHong Zhang 
318ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
3198ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
320c3d102feSKris Buschelman   } else {
321a0e1a404SHong Zhang     *newmat = B;
322c3d102feSKris Buschelman   }
323a0e1a404SHong Zhang 
324a0e1a404SHong Zhang   PetscFunctionReturn(0);
325a0e1a404SHong Zhang }
326a0e1a404SHong Zhang EXTERN_C_END
327