xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
259557b74SHong Zhang 
37c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
47c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
57c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
659557b74SHong Zhang 
759557b74SHong Zhang EXTERN_C_BEGIN
859557b74SHong Zhang #undef __FUNCT__
94124e67cSShri Abhyankar #define __FUNCT__ "MatConvert_SeqSBAIJ_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 */
2274ed9c26SBarry Smith   ierr = PetscMalloc2(m,PetscInt,&rowlengths,m+1,PetscInt,&rowstart);CHKERRQ(ierr);
23a7a3a9ebSHong Zhang 
24a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
254e5e7fe4SHong Zhang   aj = a->j;
26a7a3a9ebSHong Zhang   k = 0;
27a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
284e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
29a7a3a9ebSHong Zhang     aj++; /* skip diagonal */
30a7a3a9ebSHong Zhang     for (j=1; j<nz; j++) { /* no. of lower triangular blocks */
31a7a3a9ebSHong Zhang       rowlengths[(*aj)*bs]++; aj++;
32a7a3a9ebSHong Zhang     }
33a7a3a9ebSHong Zhang     rowlengths[k] += nz;   /* no. of upper triangular blocks */
34a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
35a7a3a9ebSHong Zhang     for (j=1; j<bs; j++) {
36a7a3a9ebSHong Zhang       rowlengths[k+j] = rowlengths[k];
37a7a3a9ebSHong Zhang     }
38a7a3a9ebSHong Zhang     k += bs;
39a7a3a9ebSHong Zhang     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
404e5e7fe4SHong Zhang   }
414e5e7fe4SHong Zhang 
427adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
43f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
44f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
45f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
464e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
47d0f46423SBarry Smith   B->rmap->bs = A->rmap->bs;
484e5e7fe4SHong Zhang 
494e5e7fe4SHong Zhang   b  = (Mat_SeqAIJ*)(B->data);
504e5e7fe4SHong Zhang   bi = b->i;
514e5e7fe4SHong Zhang   bj = b->j;
524e5e7fe4SHong Zhang   bv = b->a;
534e5e7fe4SHong Zhang 
544e5e7fe4SHong Zhang   /* set b->i */
55a7a3a9ebSHong Zhang   bi[0] = 0; rowstart[0] = 0;
56a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++){
57a7a3a9ebSHong Zhang     for (j=0; j<bs; j++){
58a7a3a9ebSHong Zhang       b->ilen[i*bs+j]  = rowlengths[i*bs];
59a7a3a9ebSHong Zhang       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
604e5e7fe4SHong Zhang     }
61a7a3a9ebSHong Zhang     bi[i+1]     = bi[i] + rowlengths[i*bs]/bs;
62a7a3a9ebSHong Zhang   }
63*e32f2f54SBarry 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);
644e5e7fe4SHong Zhang 
654e5e7fe4SHong Zhang   /* set b->j and b->a */
664e5e7fe4SHong Zhang   aj = a->j; av = a->a;
67a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
68a7a3a9ebSHong Zhang     /* diagonal block */
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;
78a7a3a9ebSHong Zhang 
794e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i] -1;
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;
86a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
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 
106ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1078ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
108c3d102feSKris Buschelman   } else {
1094e5e7fe4SHong Zhang     *newmat = B;
110c3d102feSKris Buschelman   }
1114e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1124e5e7fe4SHong Zhang }
113be1d678aSKris Buschelman EXTERN_C_END
114be1d678aSKris Buschelman 
115be1d678aSKris Buschelman EXTERN_C_BEGIN
1164e5e7fe4SHong Zhang #undef __FUNCT__
11759557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
1188cf70c4bSSatish Balay PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) {
119676c34cdSKris Buschelman   Mat            B;
12059557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
121861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
122dfbe8321SBarry Smith   PetscErrorCode ierr;
123d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
124dd6ea824SBarry Smith   MatScalar      *av,*bv;
12559557b74SHong Zhang 
12659557b74SHong Zhang   PetscFunctionBegin;
127*e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
12859557b74SHong Zhang 
12913f74950SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13059557b74SHong Zhang   for (i=0; i<m; i++) {
13159557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13259557b74SHong Zhang   }
1337adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
134f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
135f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
136ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
13759557b74SHong Zhang 
1384e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
13959557b74SHong Zhang 
140676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
141861ba921SHong Zhang   bi = b->i;
142861ba921SHong Zhang   bj = b->j;
143861ba921SHong Zhang   bv = b->a;
144861ba921SHong Zhang 
145861ba921SHong Zhang   bi[0] = 0;
14659557b74SHong Zhang   for (i=0; i<m; i++) {
14759557b74SHong Zhang     aj = a->j + a->diag[i];
14859557b74SHong Zhang     av = a->a + a->diag[i];
149861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++){
150861ba921SHong Zhang       *bj = *aj; bj++; aj++;
151861ba921SHong Zhang       *bv = *av; bv++; av++;
152861ba921SHong Zhang     }
153861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
154861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
15559557b74SHong Zhang   }
15659557b74SHong Zhang 
15759557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
158676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160eeffb40dSHong Zhang   if (A->hermitian){
161eeffb40dSHong Zhang     ierr = MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
162eeffb40dSHong Zhang   }
163676c34cdSKris Buschelman 
164ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1658ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
166c3d102feSKris Buschelman   } else {
167676c34cdSKris Buschelman     *newmat = B;
168c3d102feSKris Buschelman   }
16959557b74SHong Zhang   PetscFunctionReturn(0);
17059557b74SHong Zhang }
17159557b74SHong Zhang EXTERN_C_END
17259557b74SHong Zhang 
173a0e1a404SHong Zhang EXTERN_C_BEGIN
174a0e1a404SHong Zhang #undef __FUNCT__
175a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ"
176f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
177a0e1a404SHong Zhang {
178a0e1a404SHong Zhang   Mat            B;
179a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
180a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
181dfbe8321SBarry Smith   PetscErrorCode ierr;
182d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
183d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs;
184dd6ea824SBarry Smith   MatScalar      *av,*bv;
185a0e1a404SHong Zhang 
186a0e1a404SHong Zhang   PetscFunctionBegin;
187a0e1a404SHong Zhang   /* compute browlengths of newmat */
18874ed9c26SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);CHKERRQ(ierr);
189a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
190a0e1a404SHong Zhang   aj = a->j;
191a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
192a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
193a0e1a404SHong Zhang     aj++; /* skip diagonal */
194a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
195a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
196a0e1a404SHong Zhang     }
197a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
198a0e1a404SHong Zhang   }
199a0e1a404SHong Zhang 
2007adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
201f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
202f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
203f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
2044e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
205a0e1a404SHong Zhang 
206a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
207a0e1a404SHong Zhang   bi = b->i;
208a0e1a404SHong Zhang   bj = b->j;
209a0e1a404SHong Zhang   bv = b->a;
210a0e1a404SHong Zhang 
211a0e1a404SHong Zhang   /* set b->i */
212a0e1a404SHong Zhang   bi[0] = 0;
213a0e1a404SHong Zhang   for (i=0; i<mbs; i++){
214a0e1a404SHong Zhang     b->ilen[i]    = browlengths[i];
215a0e1a404SHong Zhang     bi[i+1]       = bi[i] + browlengths[i];
216a0e1a404SHong Zhang     browstart[i]  = bi[i];
217a0e1a404SHong Zhang   }
218*e32f2f54SBarry 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);
219a0e1a404SHong Zhang 
220a0e1a404SHong Zhang   /* set b->j and b->a */
221a0e1a404SHong Zhang   aj = a->j; av = a->a;
222a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
223a0e1a404SHong Zhang     /* diagonal block */
224a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
225a0e1a404SHong Zhang     itmp = bs2*browstart[i];
226a0e1a404SHong Zhang     for (k=0; k<bs2; k++){
227a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
228a0e1a404SHong Zhang     }
229a0e1a404SHong Zhang     browstart[i]++;
230a0e1a404SHong Zhang 
231a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
232a0e1a404SHong Zhang     while (nz--){
233a0e1a404SHong Zhang       /* lower triangular blocks */
234a0e1a404SHong Zhang       *(bj + browstart[*aj]) = i;
235a0e1a404SHong Zhang       itmp = bs2*browstart[*aj];
236a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
237a0e1a404SHong Zhang         *(bv + itmp + k) = *(av + k);
238a0e1a404SHong Zhang       }
239a0e1a404SHong Zhang       browstart[*aj]++;
240a0e1a404SHong Zhang 
241a0e1a404SHong Zhang       /* upper triangular blocks */
242a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
243a0e1a404SHong Zhang       itmp = bs2*browstart[i];
244a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
245a0e1a404SHong Zhang         *(bv + itmp + k) = *av; av++;
246a0e1a404SHong Zhang       }
247a0e1a404SHong Zhang       browstart[i]++;
248a0e1a404SHong Zhang     }
249a0e1a404SHong Zhang   }
25074ed9c26SBarry Smith   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
251a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253a0e1a404SHong Zhang 
254ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
2558ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
256c3d102feSKris Buschelman   } else {
257a0e1a404SHong Zhang     *newmat = B;
258c3d102feSKris Buschelman   }
259a0e1a404SHong Zhang   PetscFunctionReturn(0);
260a0e1a404SHong Zhang }
261be1d678aSKris Buschelman EXTERN_C_END
262be1d678aSKris Buschelman 
263be1d678aSKris Buschelman EXTERN_C_BEGIN
264a0e1a404SHong Zhang #undef __FUNCT__
265a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
266f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
267a0e1a404SHong Zhang {
268a0e1a404SHong Zhang   Mat            B;
269a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
270a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
271dfbe8321SBarry Smith   PetscErrorCode ierr;
272d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
273d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
274dd6ea824SBarry Smith   MatScalar      *av,*bv;
2752af78befSBarry Smith   PetscTruth     flg;
276a0e1a404SHong Zhang 
277a0e1a404SHong Zhang   PetscFunctionBegin;
278*e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2792af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
280*e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
281a0e1a404SHong Zhang 
28213f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
283a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
284a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
285a0e1a404SHong Zhang   }
286a0e1a404SHong Zhang 
2877adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
288f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
289f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
290ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr);
2914e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
292a0e1a404SHong Zhang 
293a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
294a0e1a404SHong Zhang   bi = b->i;
295a0e1a404SHong Zhang   bj = b->j;
296a0e1a404SHong Zhang   bv = b->a;
297a0e1a404SHong Zhang 
298a0e1a404SHong Zhang   bi[0] = 0;
299a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
300a0e1a404SHong Zhang     aj = a->j + a->diag[i];
301a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
302a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++){
303a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
304a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
305a0e1a404SHong Zhang         *bv = *av; bv++; av++;
306a0e1a404SHong Zhang       }
307a0e1a404SHong Zhang     }
308a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
309a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
310a0e1a404SHong Zhang   }
311a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
312a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
314a0e1a404SHong Zhang 
315ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
3168ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
317c3d102feSKris Buschelman   } else {
318a0e1a404SHong Zhang     *newmat = B;
319c3d102feSKris Buschelman   }
320a0e1a404SHong Zhang 
321a0e1a404SHong Zhang   PetscFunctionReturn(0);
322a0e1a404SHong Zhang }
323a0e1a404SHong Zhang EXTERN_C_END
324