xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision 367daffbeca5f1514f664d9be3910d7b6e19736f)
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"
8cc2e6a90SBarry 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;
23a7a3a9ebSHong Zhang   k  = 0;
24a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
254e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
2601be0148SBarry Smith     if (nz) {
2701be0148SBarry Smith       rowlengths[k] += nz;   /* no. of upper triangular blocks */
2801be0148SBarry Smith       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
2901be0148SBarry Smith       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
30a7a3a9ebSHong Zhang         rowlengths[(*aj)*bs]++; aj++;
31a7a3a9ebSHong Zhang       }
3201be0148SBarry Smith     }
33a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
34a7a3a9ebSHong Zhang     for (j=1; j<bs; j++) {
35a7a3a9ebSHong Zhang       rowlengths[k+j] = rowlengths[k];
36a7a3a9ebSHong Zhang     }
37a7a3a9ebSHong Zhang     k += bs;
38a7a3a9ebSHong Zhang     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
394e5e7fe4SHong Zhang   }
404e5e7fe4SHong Zhang 
41ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
42f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
436d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
44f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
454e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
4626fbe8dcSKarl Rupp 
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   }
6301be0148SBarry 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);
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++) {
6801be0148SBarry Smith     nz = ai[i+1] - ai[i];
69a7a3a9ebSHong Zhang     /* diagonal block */
7001be0148SBarry Smith     if (nz && *aj == i) {
7101be0148SBarry Smith       nz--;
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;
8101be0148SBarry Smith     }
82a7a3a9ebSHong Zhang 
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;
893bededecSBarry Smith           *(bv + rowstart[itmp]) = *(av+j*bs+k);
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   }
10574ed9c26SBarry Smith   ierr = PetscFree2(rowlengths,rowstart);CHKERRQ(ierr);
1064e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1074e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1084e5e7fe4SHong Zhang 
109511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
11028be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
111c3d102feSKris Buschelman   } else {
1124e5e7fe4SHong Zhang     *newmat = B;
113c3d102feSKris Buschelman   }
1144e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1154e5e7fe4SHong Zhang }
116be1d678aSKris Buschelman 
1174e5e7fe4SHong Zhang #undef __FUNCT__
11859557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
119cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1200adebc6cSBarry Smith {
121676c34cdSKris Buschelman   Mat            B;
12259557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
123861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
124dfbe8321SBarry Smith   PetscErrorCode ierr;
125d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
126dd6ea824SBarry Smith   MatScalar      *av,*bv;
12759557b74SHong Zhang 
12859557b74SHong Zhang   PetscFunctionBegin;
129ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
130e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
13159557b74SHong Zhang 
132785e854fSJed Brown   ierr = PetscMalloc1(m,&rowlengths);CHKERRQ(ierr);
13359557b74SHong Zhang   for (i=0; i<m; i++) {
13459557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13559557b74SHong Zhang   }
136ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
137f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1386d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
139*367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr);
14059557b74SHong Zhang 
1414e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
14259557b74SHong Zhang 
143676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
144861ba921SHong Zhang   bi = b->i;
145861ba921SHong Zhang   bj = b->j;
146861ba921SHong Zhang   bv = b->a;
147861ba921SHong Zhang 
148861ba921SHong Zhang   bi[0] = 0;
14959557b74SHong Zhang   for (i=0; i<m; i++) {
15059557b74SHong Zhang     aj = a->j + a->diag[i];
15159557b74SHong Zhang     av = a->a + a->diag[i];
152861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++) {
153861ba921SHong Zhang       *bj = *aj; bj++; aj++;
154861ba921SHong Zhang       *bv = *av; bv++; av++;
155861ba921SHong Zhang     }
156861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
157861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
15859557b74SHong Zhang   }
15959557b74SHong Zhang 
16059557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
161676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163676c34cdSKris Buschelman 
164511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16528be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
166c3d102feSKris Buschelman   } else {
167676c34cdSKris Buschelman     *newmat = B;
168c3d102feSKris Buschelman   }
16959557b74SHong Zhang   PetscFunctionReturn(0);
17059557b74SHong Zhang }
17159557b74SHong Zhang 
172a0e1a404SHong Zhang #undef __FUNCT__
1737af9e4fdSHong Zhang #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqBAIJ"
174cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode 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;
18174ee4d9fSHong Zhang   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
182dd6ea824SBarry Smith   MatScalar      *av,*bv;
183a0e1a404SHong Zhang 
184a0e1a404SHong Zhang   PetscFunctionBegin;
185a0e1a404SHong Zhang   /* compute browlengths of newmat */
186dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&browlengths,mbs,&browstart);CHKERRQ(ierr);
187a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
188a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
189a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
190a0e1a404SHong Zhang     aj++; /* skip diagonal */
191a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
192a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
193a0e1a404SHong Zhang     }
194a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
195a0e1a404SHong Zhang   }
196a0e1a404SHong Zhang 
197ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
198f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1996d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
200f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
2014e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
202a0e1a404SHong Zhang 
203a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
204a0e1a404SHong Zhang   bi = b->i;
205a0e1a404SHong Zhang   bj = b->j;
206a0e1a404SHong Zhang   bv = b->a;
207a0e1a404SHong Zhang 
208a0e1a404SHong Zhang   /* set b->i */
209a0e1a404SHong Zhang   bi[0] = 0;
210a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
211a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
212a0e1a404SHong Zhang     bi[i+1]      = bi[i] + browlengths[i];
213a0e1a404SHong Zhang     browstart[i] = bi[i];
214a0e1a404SHong Zhang   }
215e32f2f54SBarry 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);
216a0e1a404SHong Zhang 
217a0e1a404SHong Zhang   /* set b->j and b->a */
218a0e1a404SHong Zhang   aj = a->j; av = a->a;
219a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
220a0e1a404SHong Zhang     /* diagonal block */
221a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
22226fbe8dcSKarl Rupp 
223a0e1a404SHong Zhang     itmp = bs2*browstart[i];
224a0e1a404SHong Zhang     for (k=0; k<bs2; k++) {
225a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
226a0e1a404SHong Zhang     }
227a0e1a404SHong Zhang     browstart[i]++;
228a0e1a404SHong Zhang 
229a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
230a0e1a404SHong Zhang     while (nz--) {
23174ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
23274ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
23326fbe8dcSKarl Rupp 
23474ee4d9fSHong Zhang       itmp = bs2*browstart[*aj];  /* row index */
23574ee4d9fSHong Zhang       for (col=0; col<bs; col++) {
23674ee4d9fSHong Zhang         k = col;
23774ee4d9fSHong Zhang         for (row=0; row<bs; row++) {
23874ee4d9fSHong Zhang           bv[itmp + col*bs+row] = av[k]; k+=bs;
23974ee4d9fSHong Zhang         }
240a0e1a404SHong Zhang       }
241a0e1a404SHong Zhang       browstart[*aj]++;
242a0e1a404SHong Zhang 
243a0e1a404SHong Zhang       /* upper triangular blocks */
244a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
24526fbe8dcSKarl Rupp 
246a0e1a404SHong Zhang       itmp = bs2*browstart[i];
247a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
24874ee4d9fSHong Zhang         bv[itmp + k] = av[k];
249a0e1a404SHong Zhang       }
25074ee4d9fSHong Zhang       av += bs2;
251a0e1a404SHong Zhang       browstart[i]++;
252a0e1a404SHong Zhang     }
253a0e1a404SHong Zhang   }
25474ed9c26SBarry Smith   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
255a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257a0e1a404SHong Zhang 
258511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
25928be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
260c3d102feSKris Buschelman   } else {
261a0e1a404SHong Zhang     *newmat = B;
262c3d102feSKris Buschelman   }
263a0e1a404SHong Zhang   PetscFunctionReturn(0);
264a0e1a404SHong Zhang }
265be1d678aSKris Buschelman 
266a0e1a404SHong Zhang #undef __FUNCT__
267a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
268cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
269a0e1a404SHong Zhang {
270a0e1a404SHong Zhang   Mat            B;
271a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
272a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
273dfbe8321SBarry Smith   PetscErrorCode ierr;
274d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
275d0f46423SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
276dd6ea824SBarry Smith   MatScalar      *av,*bv;
277ace3abfcSBarry Smith   PetscBool      flg;
278a0e1a404SHong Zhang 
279a0e1a404SHong Zhang   PetscFunctionBegin;
280ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
281e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2822af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
283e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
284a0e1a404SHong Zhang 
285785e854fSJed Brown   ierr = PetscMalloc1(mbs,&browlengths);CHKERRQ(ierr);
286a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
287a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
288a0e1a404SHong Zhang   }
289a0e1a404SHong Zhang 
290ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
291f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2926d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
293*367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
2944e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);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 
318511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
31928be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
320c3d102feSKris Buschelman   } else {
321a0e1a404SHong Zhang     *newmat = B;
322c3d102feSKris Buschelman   }
323a0e1a404SHong Zhang   PetscFunctionReturn(0);
324a0e1a404SHong Zhang }
325