xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision cc2e6a90c05b27ffec69cb207fe793d447f14420)
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 
659557b74SHong Zhang #undef __FUNCT__
74124e67cSShri Abhyankar #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqAIJ"
8*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
94e5e7fe4SHong Zhang {
104e5e7fe4SHong Zhang   Mat            B;
114e5e7fe4SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
124e5e7fe4SHong Zhang   Mat_SeqAIJ     *b;
13dfbe8321SBarry Smith   PetscErrorCode ierr;
14d0f46423SBarry 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;
1501be0148SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0;
16dd6ea824SBarry Smith   MatScalar      *av,*bv;
174e5e7fe4SHong Zhang 
184e5e7fe4SHong Zhang   PetscFunctionBegin;
194e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
20dcca6d9dSJed Brown   ierr = PetscMalloc2(m,&rowlengths,m+1,&rowstart);CHKERRQ(ierr);
21a7a3a9ebSHong Zhang 
22a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
234e5e7fe4SHong Zhang   aj = a->j;
24a7a3a9ebSHong Zhang   k  = 0;
25a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
264e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
2701be0148SBarry Smith     if (nz) {
2801be0148SBarry Smith       rowlengths[k] += nz;   /* no. of upper triangular blocks */
2901be0148SBarry Smith       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
3001be0148SBarry Smith       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
31a7a3a9ebSHong Zhang         rowlengths[(*aj)*bs]++; aj++;
32a7a3a9ebSHong Zhang       }
3301be0148SBarry Smith     }
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 
42ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
43f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
446d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
45f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
464e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
4726fbe8dcSKarl Rupp 
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   }
6401be0148SBarry 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);
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++) {
6901be0148SBarry Smith     nz = ai[i+1] - ai[i];
70a7a3a9ebSHong Zhang     /* diagonal block */
7101be0148SBarry Smith     if (nz && *aj == i) {
7201be0148SBarry Smith       nz--;
73a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row i*bs+j */
74a7a3a9ebSHong Zhang         itmp = i*bs+j;
75a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col i*bs+k */
76a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
77a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
78a7a3a9ebSHong Zhang           rowstart[itmp]++;
79a7a3a9ebSHong Zhang         }
80a7a3a9ebSHong Zhang       }
81a7a3a9ebSHong Zhang       aj++; av += bs2;
8201be0148SBarry Smith     }
83a7a3a9ebSHong Zhang 
844e5e7fe4SHong Zhang     while (nz--) {
85a7a3a9ebSHong Zhang       /* lower triangular blocks */
86a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row (*aj)*bs+j */
87a7a3a9ebSHong Zhang         itmp = (*aj)*bs+j;
88a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col i*bs+k */
89a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = i*bs+k;
903bededecSBarry Smith           *(bv + rowstart[itmp]) = *(av+j*bs+k);
91a7a3a9ebSHong Zhang           rowstart[itmp]++;
92a7a3a9ebSHong Zhang         }
93a7a3a9ebSHong Zhang       }
94a7a3a9ebSHong Zhang       /* upper triangular blocks */
95a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row i*bs+j */
96a7a3a9ebSHong Zhang         itmp = i*bs+j;
97a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col (*aj)*bs+k */
98a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
99a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
100a7a3a9ebSHong Zhang           rowstart[itmp]++;
101a7a3a9ebSHong Zhang         }
102a7a3a9ebSHong Zhang       }
103a7a3a9ebSHong Zhang       aj++; av += bs2;
1044e5e7fe4SHong Zhang     }
1054e5e7fe4SHong Zhang   }
10674ed9c26SBarry Smith   ierr = PetscFree2(rowlengths,rowstart);CHKERRQ(ierr);
1074e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1084e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1094e5e7fe4SHong Zhang 
110511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
11128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
112c3d102feSKris Buschelman   } else {
1134e5e7fe4SHong Zhang     *newmat = B;
114c3d102feSKris Buschelman   }
1154e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1164e5e7fe4SHong Zhang }
117be1d678aSKris Buschelman 
1184e5e7fe4SHong Zhang #undef __FUNCT__
11959557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
120*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1210adebc6cSBarry Smith {
122676c34cdSKris Buschelman   Mat            B;
12359557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
124861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
125dfbe8321SBarry Smith   PetscErrorCode ierr;
126d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
127dd6ea824SBarry Smith   MatScalar      *av,*bv;
12859557b74SHong Zhang 
12959557b74SHong Zhang   PetscFunctionBegin;
130ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
131e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
13259557b74SHong Zhang 
133785e854fSJed Brown   ierr = PetscMalloc1(m,&rowlengths);CHKERRQ(ierr);
13459557b74SHong Zhang   for (i=0; i<m; i++) {
13559557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13659557b74SHong Zhang   }
137ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
138f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1396d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
140ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
14159557b74SHong Zhang 
1424e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
14359557b74SHong Zhang 
144676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
145861ba921SHong Zhang   bi = b->i;
146861ba921SHong Zhang   bj = b->j;
147861ba921SHong Zhang   bv = b->a;
148861ba921SHong Zhang 
149861ba921SHong Zhang   bi[0] = 0;
15059557b74SHong Zhang   for (i=0; i<m; i++) {
15159557b74SHong Zhang     aj = a->j + a->diag[i];
15259557b74SHong Zhang     av = a->a + a->diag[i];
153861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++) {
154861ba921SHong Zhang       *bj = *aj; bj++; aj++;
155861ba921SHong Zhang       *bv = *av; bv++; av++;
156861ba921SHong Zhang     }
157861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
158861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
15959557b74SHong Zhang   }
16059557b74SHong Zhang 
16159557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
162676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164676c34cdSKris Buschelman 
165511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16628be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
167c3d102feSKris Buschelman   } else {
168676c34cdSKris Buschelman     *newmat = B;
169c3d102feSKris Buschelman   }
17059557b74SHong Zhang   PetscFunctionReturn(0);
17159557b74SHong Zhang }
17259557b74SHong Zhang 
173a0e1a404SHong Zhang #undef __FUNCT__
1747af9e4fdSHong Zhang #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqBAIJ"
175*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
176a0e1a404SHong Zhang {
177a0e1a404SHong Zhang   Mat            B;
178a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
179a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
180dfbe8321SBarry Smith   PetscErrorCode ierr;
181d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
18274ee4d9fSHong Zhang   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
183dd6ea824SBarry Smith   MatScalar      *av,*bv;
184a0e1a404SHong Zhang 
185a0e1a404SHong Zhang   PetscFunctionBegin;
186a0e1a404SHong Zhang   /* compute browlengths of newmat */
187dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&browlengths,mbs,&browstart);CHKERRQ(ierr);
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 
199ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
200f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2016d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQBAIJ);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   }
217e32f2f54SBarry 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);
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++;
22426fbe8dcSKarl Rupp 
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--) {
23374ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
23474ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
23526fbe8dcSKarl Rupp 
23674ee4d9fSHong Zhang       itmp = bs2*browstart[*aj];  /* row index */
23774ee4d9fSHong Zhang       for (col=0; col<bs; col++) {
23874ee4d9fSHong Zhang         k = col;
23974ee4d9fSHong Zhang         for (row=0; row<bs; row++) {
24074ee4d9fSHong Zhang           bv[itmp + col*bs+row] = av[k]; k+=bs;
24174ee4d9fSHong Zhang         }
242a0e1a404SHong Zhang       }
243a0e1a404SHong Zhang       browstart[*aj]++;
244a0e1a404SHong Zhang 
245a0e1a404SHong Zhang       /* upper triangular blocks */
246a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
24726fbe8dcSKarl Rupp 
248a0e1a404SHong Zhang       itmp = bs2*browstart[i];
249a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
25074ee4d9fSHong Zhang         bv[itmp + k] = av[k];
251a0e1a404SHong Zhang       }
25274ee4d9fSHong Zhang       av += bs2;
253a0e1a404SHong Zhang       browstart[i]++;
254a0e1a404SHong Zhang     }
255a0e1a404SHong Zhang   }
25674ed9c26SBarry Smith   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
257a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259a0e1a404SHong Zhang 
260511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
26128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
262c3d102feSKris Buschelman   } else {
263a0e1a404SHong Zhang     *newmat = B;
264c3d102feSKris Buschelman   }
265a0e1a404SHong Zhang   PetscFunctionReturn(0);
266a0e1a404SHong Zhang }
267be1d678aSKris Buschelman 
268a0e1a404SHong Zhang #undef __FUNCT__
269a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
270*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, 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;
276d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
277d0f46423SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
278dd6ea824SBarry Smith   MatScalar      *av,*bv;
279ace3abfcSBarry Smith   PetscBool      flg;
280a0e1a404SHong Zhang 
281a0e1a404SHong Zhang   PetscFunctionBegin;
282ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
283e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2842af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
285e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
286a0e1a404SHong Zhang 
287785e854fSJed Brown   ierr = PetscMalloc1(mbs,&browlengths);CHKERRQ(ierr);
288a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
289a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
290a0e1a404SHong Zhang   }
291a0e1a404SHong Zhang 
292ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
293f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2946d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
295ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr);
2964e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
297a0e1a404SHong Zhang 
298a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
299a0e1a404SHong Zhang   bi = b->i;
300a0e1a404SHong Zhang   bj = b->j;
301a0e1a404SHong Zhang   bv = b->a;
302a0e1a404SHong Zhang 
303a0e1a404SHong Zhang   bi[0] = 0;
304a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
305a0e1a404SHong Zhang     aj = a->j + a->diag[i];
306a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
307a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++) {
308a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
309a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
310a0e1a404SHong Zhang         *bv = *av; bv++; av++;
311a0e1a404SHong Zhang       }
312a0e1a404SHong Zhang     }
313a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
314a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
315a0e1a404SHong Zhang   }
316a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
317a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
318a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
319a0e1a404SHong Zhang 
320511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
32128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
322c3d102feSKris Buschelman   } else {
323a0e1a404SHong Zhang     *newmat = B;
324c3d102feSKris Buschelman   }
325a0e1a404SHong Zhang   PetscFunctionReturn(0);
326a0e1a404SHong Zhang }
327