xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision eb1ec7c1db909b5b8e5e7107245f4f2eaa6849e7)
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 
6cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
74e5e7fe4SHong Zhang {
84e5e7fe4SHong Zhang   Mat            B;
94e5e7fe4SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
104e5e7fe4SHong Zhang   Mat_SeqAIJ     *b;
11dfbe8321SBarry Smith   PetscErrorCode ierr;
12d0f46423SBarry 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;
1301be0148SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0;
14dd6ea824SBarry Smith   MatScalar      *av,*bv;
15*eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
16*eb1ec7c1SStefano Zampini   const int      aconj = A->hermitian ? 1 : 0;
17*eb1ec7c1SStefano Zampini #else
18*eb1ec7c1SStefano Zampini   const int      aconj = 0;
19*eb1ec7c1SStefano Zampini #endif
204e5e7fe4SHong Zhang 
214e5e7fe4SHong Zhang   PetscFunctionBegin;
224e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
23dcca6d9dSJed Brown   ierr = PetscMalloc2(m,&rowlengths,m+1,&rowstart);CHKERRQ(ierr);
24a7a3a9ebSHong Zhang 
25a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
26a7a3a9ebSHong Zhang   k  = 0;
27a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
284e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
2901be0148SBarry Smith     if (nz) {
3001be0148SBarry Smith       rowlengths[k] += nz;   /* no. of upper triangular blocks */
3101be0148SBarry Smith       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
3201be0148SBarry Smith       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
33a7a3a9ebSHong Zhang         rowlengths[(*aj)*bs]++; aj++;
34a7a3a9ebSHong Zhang       }
3501be0148SBarry Smith     }
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;
414e5e7fe4SHong Zhang   }
424e5e7fe4SHong Zhang 
43bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
44ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
45f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
466d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
47f204ca49SKris Buschelman     ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
48bd019fc1SStefano Zampini     ierr = MatSetBlockSize(B,A->rmap->bs);CHKERRQ(ierr);
49bd019fc1SStefano Zampini   } else B = *newmat;
504e5e7fe4SHong Zhang 
514e5e7fe4SHong Zhang   b  = (Mat_SeqAIJ*)(B->data);
524e5e7fe4SHong Zhang   bi = b->i;
534e5e7fe4SHong Zhang   bj = b->j;
544e5e7fe4SHong Zhang   bv = b->a;
554e5e7fe4SHong Zhang 
564e5e7fe4SHong Zhang   /* set b->i */
57a7a3a9ebSHong Zhang   bi[0] = 0; rowstart[0] = 0;
58a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
59a7a3a9ebSHong Zhang     for (j=0; j<bs; j++) {
60a7a3a9ebSHong Zhang       b->ilen[i*bs+j]    = rowlengths[i*bs];
61a7a3a9ebSHong Zhang       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
624e5e7fe4SHong Zhang     }
63a7a3a9ebSHong Zhang     bi[i+1] = bi[i] + rowlengths[i*bs]/bs;
64a7a3a9ebSHong Zhang   }
6501be0148SBarry 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);
664e5e7fe4SHong Zhang 
674e5e7fe4SHong Zhang   /* set b->j and b->a */
684e5e7fe4SHong Zhang   aj = a->j; av = a->a;
69a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
7001be0148SBarry Smith     nz = ai[i+1] - ai[i];
71a7a3a9ebSHong Zhang     /* diagonal block */
7201be0148SBarry Smith     if (nz && *aj == i) {
7301be0148SBarry Smith       nz--;
74a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row i*bs+j */
75a7a3a9ebSHong Zhang         itmp = i*bs+j;
76a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col i*bs+k */
77a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
78a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
79a7a3a9ebSHong Zhang           rowstart[itmp]++;
80a7a3a9ebSHong Zhang         }
81a7a3a9ebSHong Zhang       }
82a7a3a9ebSHong Zhang       aj++; av += bs2;
8301be0148SBarry Smith     }
84a7a3a9ebSHong Zhang 
854e5e7fe4SHong Zhang     while (nz--) {
86a7a3a9ebSHong Zhang       /* lower triangular blocks */
87a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row (*aj)*bs+j */
88a7a3a9ebSHong Zhang         itmp = (*aj)*bs+j;
89a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col i*bs+k */
90a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = i*bs+k;
91*eb1ec7c1SStefano Zampini           *(bv + rowstart[itmp]) = aconj ? PetscConj(*(av+j*bs+k)) : *(av+j*bs+k);
92a7a3a9ebSHong Zhang           rowstart[itmp]++;
93a7a3a9ebSHong Zhang         }
94a7a3a9ebSHong Zhang       }
95a7a3a9ebSHong Zhang       /* upper triangular blocks */
96a7a3a9ebSHong Zhang       for (j=0; j<bs; j++) {   /* row i*bs+j */
97a7a3a9ebSHong Zhang         itmp = i*bs+j;
98a7a3a9ebSHong Zhang         for (k=0; k<bs; k++) { /* col (*aj)*bs+k */
99a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj)*bs+k;
100a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
101a7a3a9ebSHong Zhang           rowstart[itmp]++;
102a7a3a9ebSHong Zhang         }
103a7a3a9ebSHong Zhang       }
104a7a3a9ebSHong Zhang       aj++; av += bs2;
1054e5e7fe4SHong Zhang     }
1064e5e7fe4SHong Zhang   }
10774ed9c26SBarry Smith   ierr = PetscFree2(rowlengths,rowstart);CHKERRQ(ierr);
1084e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1094e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1104e5e7fe4SHong Zhang 
111511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
11228be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
113c3d102feSKris Buschelman   } else {
1144e5e7fe4SHong Zhang     *newmat = B;
115c3d102feSKris Buschelman   }
1164e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1174e5e7fe4SHong Zhang }
118be1d678aSKris Buschelman 
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;
125ae8d29abSPierre Jolivet   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths,bs=PetscAbs(A->rmap->bs);
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 
132ae8d29abSPierre Jolivet   ierr = PetscMalloc1(m/bs,&rowlengths);CHKERRQ(ierr);
133ae8d29abSPierre Jolivet   for (i=0; i<m/bs; i++) {
134ae8d29abSPierre Jolivet     rowlengths[i] = (ai[i*bs+1] - a->diag[i*bs])/bs;
13559557b74SHong Zhang   }
136bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
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);
140ae8d29abSPierre Jolivet     ierr = MatSeqSBAIJSetPreallocation(B,bs,0,rowlengths);CHKERRQ(ierr);
141bd019fc1SStefano Zampini   } else B = *newmat;
14259557b74SHong Zhang 
143ae8d29abSPierre Jolivet   if (bs == 1) {
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     }
160676c34cdSKris Buschelman     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161676c34cdSKris Buschelman     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162ae8d29abSPierre Jolivet   } else {
163ae8d29abSPierre Jolivet     ierr = MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
164ae8d29abSPierre Jolivet     /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
165ae8d29abSPierre Jolivet     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
166ae8d29abSPierre Jolivet     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
167ae8d29abSPierre Jolivet     ierr = MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&B);CHKERRQ(ierr);
168ae8d29abSPierre Jolivet   }
169ae8d29abSPierre Jolivet   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
170511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
17128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
172ae8d29abSPierre Jolivet   } else *newmat = B;
17359557b74SHong Zhang   PetscFunctionReturn(0);
17459557b74SHong Zhang }
17559557b74SHong Zhang 
176cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode 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;
18374ee4d9fSHong Zhang   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
184dd6ea824SBarry Smith   MatScalar      *av,*bv;
185*eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
186*eb1ec7c1SStefano Zampini   const int      aconj = A->hermitian ? 1 : 0;
187*eb1ec7c1SStefano Zampini #else
188*eb1ec7c1SStefano Zampini   const int      aconj = 0;
189*eb1ec7c1SStefano Zampini #endif
190a0e1a404SHong Zhang 
191a0e1a404SHong Zhang   PetscFunctionBegin;
192a0e1a404SHong Zhang   /* compute browlengths of newmat */
193dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&browlengths,mbs,&browstart);CHKERRQ(ierr);
194a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
195a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
196a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
197a0e1a404SHong Zhang     aj++; /* skip diagonal */
198a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
199a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
200a0e1a404SHong Zhang     }
201a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
202a0e1a404SHong Zhang   }
203a0e1a404SHong Zhang 
204bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
205ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
206f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2076d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
208f204ca49SKris Buschelman     ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
209bd019fc1SStefano Zampini   } else B = *newmat;
210a0e1a404SHong Zhang 
211a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
212a0e1a404SHong Zhang   bi = b->i;
213a0e1a404SHong Zhang   bj = b->j;
214a0e1a404SHong Zhang   bv = b->a;
215a0e1a404SHong Zhang 
216a0e1a404SHong Zhang   /* set b->i */
217a0e1a404SHong Zhang   bi[0] = 0;
218a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
219a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
220a0e1a404SHong Zhang     bi[i+1]      = bi[i] + browlengths[i];
221a0e1a404SHong Zhang     browstart[i] = bi[i];
222a0e1a404SHong Zhang   }
223e32f2f54SBarry 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);
224a0e1a404SHong Zhang 
225a0e1a404SHong Zhang   /* set b->j and b->a */
226a0e1a404SHong Zhang   aj = a->j; av = a->a;
227a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
228a0e1a404SHong Zhang     /* diagonal block */
229a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
23026fbe8dcSKarl Rupp 
231a0e1a404SHong Zhang     itmp = bs2*browstart[i];
232a0e1a404SHong Zhang     for (k=0; k<bs2; k++) {
233a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
234a0e1a404SHong Zhang     }
235a0e1a404SHong Zhang     browstart[i]++;
236a0e1a404SHong Zhang 
237a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
238a0e1a404SHong Zhang     while (nz--) {
23974ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
24074ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
24126fbe8dcSKarl Rupp 
24274ee4d9fSHong Zhang       itmp = bs2*browstart[*aj];  /* row index */
24374ee4d9fSHong Zhang       for (col=0; col<bs; col++) {
24474ee4d9fSHong Zhang         k = col;
24574ee4d9fSHong Zhang         for (row=0; row<bs; row++) {
246*eb1ec7c1SStefano Zampini           bv[itmp + col*bs+row] = aconj ? PetscConj(av[k]) : av[k];
247*eb1ec7c1SStefano Zampini           k+=bs;
24874ee4d9fSHong Zhang         }
249a0e1a404SHong Zhang       }
250a0e1a404SHong Zhang       browstart[*aj]++;
251a0e1a404SHong Zhang 
252a0e1a404SHong Zhang       /* upper triangular blocks */
253a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
25426fbe8dcSKarl Rupp 
255a0e1a404SHong Zhang       itmp = bs2*browstart[i];
256a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
25774ee4d9fSHong Zhang         bv[itmp + k] = av[k];
258a0e1a404SHong Zhang       }
25974ee4d9fSHong Zhang       av += bs2;
260a0e1a404SHong Zhang       browstart[i]++;
261a0e1a404SHong Zhang     }
262a0e1a404SHong Zhang   }
26374ed9c26SBarry Smith   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
264a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
265a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
266a0e1a404SHong Zhang 
267511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
26828be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
269ae8d29abSPierre Jolivet   } else *newmat = B;
270a0e1a404SHong Zhang   PetscFunctionReturn(0);
271a0e1a404SHong Zhang }
272be1d678aSKris Buschelman 
273cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
274a0e1a404SHong Zhang {
275a0e1a404SHong Zhang   Mat            B;
276a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
277a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
280d0f46423SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
281dd6ea824SBarry Smith   MatScalar      *av,*bv;
282ace3abfcSBarry Smith   PetscBool      flg;
283a0e1a404SHong Zhang 
284a0e1a404SHong Zhang   PetscFunctionBegin;
285ce94432eSBarry Smith   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
286e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2872af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
288e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
289a0e1a404SHong Zhang 
290785e854fSJed Brown   ierr = PetscMalloc1(mbs,&browlengths);CHKERRQ(ierr);
291a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
292a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
293a0e1a404SHong Zhang   }
294a0e1a404SHong Zhang 
295bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
296ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
297f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2986d0a4a0eSHong Zhang     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
299367daffbSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
300bd019fc1SStefano Zampini   } else B = *newmat;
301a0e1a404SHong Zhang 
302a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
303a0e1a404SHong Zhang   bi = b->i;
304a0e1a404SHong Zhang   bj = b->j;
305a0e1a404SHong Zhang   bv = b->a;
306a0e1a404SHong Zhang 
307a0e1a404SHong Zhang   bi[0] = 0;
308a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
309a0e1a404SHong Zhang     aj = a->j + a->diag[i];
310a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
311a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++) {
312a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
313a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
314a0e1a404SHong Zhang         *bv = *av; bv++; av++;
315a0e1a404SHong Zhang       }
316a0e1a404SHong Zhang     }
317a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
318a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
319a0e1a404SHong Zhang   }
320a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
321a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
322a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
323a0e1a404SHong Zhang 
324511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
32528be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
326ae8d29abSPierre Jolivet   } else *newmat = B;
327a0e1a404SHong Zhang   PetscFunctionReturn(0);
328a0e1a404SHong Zhang }
329