xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision f69a0ea3504314d029ee423e0de2950ece2dab86)
1be1d678aSKris 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*f69a0ea3SMatthew 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;
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 
44*f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
45*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
46f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
47f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
48a7a3a9ebSHong Zhang   ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
494e5e7fe4SHong Zhang   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
504e5e7fe4SHong Zhang   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
51521d7252SBarry Smith   B->bs = A->bs;
524e5e7fe4SHong Zhang 
534e5e7fe4SHong Zhang   b  = (Mat_SeqAIJ*)(B->data);
544e5e7fe4SHong Zhang   bi = b->i;
554e5e7fe4SHong Zhang   bj = b->j;
564e5e7fe4SHong Zhang   bv = b->a;
574e5e7fe4SHong Zhang 
584e5e7fe4SHong Zhang   /* set b->i */
59a7a3a9ebSHong Zhang   bi[0] = 0; rowstart[0] = 0;
60a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++){
61a7a3a9ebSHong Zhang     for (j=0; j<bs; j++){
62a7a3a9ebSHong Zhang       b->ilen[i*bs+j]  = rowlengths[i*bs];
63a7a3a9ebSHong Zhang       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
644e5e7fe4SHong Zhang     }
65a7a3a9ebSHong Zhang     bi[i+1]     = bi[i] + rowlengths[i*bs]/bs;
66a7a3a9ebSHong Zhang   }
6777431f27SBarry 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);
684e5e7fe4SHong Zhang 
694e5e7fe4SHong Zhang   /* set b->j and b->a */
704e5e7fe4SHong Zhang   aj = a->j; av = a->a;
71a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
72a7a3a9ebSHong Zhang     /* diagonal block */
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;
82a7a3a9ebSHong Zhang 
834e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i] -1;
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;
90a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av+k*bs+j);
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   }
1064e5e7fe4SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1074e5e7fe4SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1084e5e7fe4SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1094e5e7fe4SHong Zhang 
110ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1118ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
112c3d102feSKris Buschelman   } else {
1134e5e7fe4SHong Zhang     *newmat = B;
114c3d102feSKris Buschelman   }
1154e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1164e5e7fe4SHong Zhang }
117be1d678aSKris Buschelman EXTERN_C_END
118be1d678aSKris Buschelman 
119be1d678aSKris Buschelman EXTERN_C_BEGIN
1204e5e7fe4SHong Zhang #undef __FUNCT__
12159557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
122*f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) {
123676c34cdSKris Buschelman   Mat            B;
12459557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
125861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
126dfbe8321SBarry Smith   PetscErrorCode ierr;
12779416396SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->M,n=A->N,i,j,*bi,*bj,*rowlengths;
128861ba921SHong Zhang   PetscScalar    *av,*bv;
12959557b74SHong Zhang 
13059557b74SHong Zhang   PetscFunctionBegin;
1312d9a3abdSHong Zhang   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
13218f87311SHong Zhang   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
13359557b74SHong Zhang 
13413f74950SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13559557b74SHong Zhang   for (i=0; i<m; i++) {
13659557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13759557b74SHong Zhang   }
138*f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
139*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
140f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
141ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
14259557b74SHong Zhang 
143676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
144676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
145676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
14659557b74SHong Zhang 
147676c34cdSKris Buschelman   b  = (Mat_SeqSBAIJ*)(B->data);
148861ba921SHong Zhang   bi = b->i;
149861ba921SHong Zhang   bj = b->j;
150861ba921SHong Zhang   bv = b->a;
151861ba921SHong Zhang 
152861ba921SHong Zhang   bi[0] = 0;
15359557b74SHong Zhang   for (i=0; i<m; i++) {
15459557b74SHong Zhang     aj = a->j + a->diag[i];
15559557b74SHong Zhang     av = a->a + a->diag[i];
156861ba921SHong Zhang     for (j=0; j<rowlengths[i]; j++){
157861ba921SHong Zhang       *bj = *aj; bj++; aj++;
158861ba921SHong Zhang       *bv = *av; bv++; av++;
159861ba921SHong Zhang     }
160861ba921SHong Zhang     bi[i+1]    = bi[i] + rowlengths[i];
161861ba921SHong Zhang     b->ilen[i] = rowlengths[i];
16259557b74SHong Zhang   }
16359557b74SHong Zhang 
16459557b74SHong Zhang   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
165676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167676c34cdSKris Buschelman 
168ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1698ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
170c3d102feSKris Buschelman   } else {
171676c34cdSKris Buschelman     *newmat = B;
172c3d102feSKris Buschelman   }
17359557b74SHong Zhang   PetscFunctionReturn(0);
17459557b74SHong Zhang }
17559557b74SHong Zhang EXTERN_C_END
17659557b74SHong Zhang 
177a0e1a404SHong Zhang EXTERN_C_BEGIN
178a0e1a404SHong Zhang #undef __FUNCT__
179a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ"
180*f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
181a0e1a404SHong Zhang {
182a0e1a404SHong Zhang   Mat            B;
183a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
184a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
185dfbe8321SBarry Smith   PetscErrorCode ierr;
18613f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
18713f74950SBarry Smith   PetscInt       bs=A->bs,bs2=bs*bs,mbs=m/bs;
188a0e1a404SHong Zhang   PetscScalar    *av,*bv;
189a0e1a404SHong Zhang 
190a0e1a404SHong Zhang   PetscFunctionBegin;
191a0e1a404SHong Zhang   /* compute browlengths of newmat */
19213f74950SBarry Smith   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
193a0e1a404SHong Zhang   browstart = browlengths + mbs;
194a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
195a0e1a404SHong Zhang   aj = a->j;
196a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
197a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
198a0e1a404SHong Zhang     aj++; /* skip diagonal */
199a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
200a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
201a0e1a404SHong Zhang     }
202a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
203a0e1a404SHong Zhang   }
204a0e1a404SHong Zhang 
205*f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
206*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
207f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
208f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
209a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
210a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
211a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
212a0e1a404SHong Zhang 
213a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
214a0e1a404SHong Zhang   bi = b->i;
215a0e1a404SHong Zhang   bj = b->j;
216a0e1a404SHong Zhang   bv = b->a;
217a0e1a404SHong Zhang 
218a0e1a404SHong Zhang   /* set b->i */
219a0e1a404SHong Zhang   bi[0] = 0;
220a0e1a404SHong Zhang   for (i=0; i<mbs; i++){
221a0e1a404SHong Zhang     b->ilen[i]    = browlengths[i];
222a0e1a404SHong Zhang     bi[i+1]       = bi[i] + browlengths[i];
223a0e1a404SHong Zhang     browstart[i]  = bi[i];
224a0e1a404SHong Zhang   }
22577431f27SBarry 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);
226a0e1a404SHong Zhang 
227a0e1a404SHong Zhang   /* set b->j and b->a */
228a0e1a404SHong Zhang   aj = a->j; av = a->a;
229a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
230a0e1a404SHong Zhang     /* diagonal block */
231a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
232a0e1a404SHong Zhang     itmp = bs2*browstart[i];
233a0e1a404SHong Zhang     for (k=0; k<bs2; k++){
234a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
235a0e1a404SHong Zhang     }
236a0e1a404SHong Zhang     browstart[i]++;
237a0e1a404SHong Zhang 
238a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
239a0e1a404SHong Zhang     while (nz--){
240a0e1a404SHong Zhang       /* lower triangular blocks */
241a0e1a404SHong Zhang       *(bj + browstart[*aj]) = i;
242a0e1a404SHong Zhang       itmp = bs2*browstart[*aj];
243a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
244a0e1a404SHong Zhang         *(bv + itmp + k) = *(av + k);
245a0e1a404SHong Zhang       }
246a0e1a404SHong Zhang       browstart[*aj]++;
247a0e1a404SHong Zhang 
248a0e1a404SHong Zhang       /* upper triangular blocks */
249a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
250a0e1a404SHong Zhang       itmp = bs2*browstart[i];
251a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
252a0e1a404SHong Zhang         *(bv + itmp + k) = *av; av++;
253a0e1a404SHong Zhang       }
254a0e1a404SHong Zhang       browstart[i]++;
255a0e1a404SHong Zhang     }
256a0e1a404SHong Zhang   }
257a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
258a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260a0e1a404SHong Zhang 
261ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
2628ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
263c3d102feSKris Buschelman   } else {
264a0e1a404SHong Zhang     *newmat = B;
265c3d102feSKris Buschelman   }
266a0e1a404SHong Zhang   PetscFunctionReturn(0);
267a0e1a404SHong Zhang }
268be1d678aSKris Buschelman EXTERN_C_END
269be1d678aSKris Buschelman 
270be1d678aSKris Buschelman EXTERN_C_BEGIN
271a0e1a404SHong Zhang #undef __FUNCT__
272a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
273*f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT 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;
27913f74950SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->m,n=A->n,i,j,k,*bi,*bj,*browlengths;
28013f74950SBarry Smith   PetscInt       bs=A->bs,bs2=bs*bs,mbs=m/bs;
281a0e1a404SHong Zhang   PetscScalar    *av,*bv;
282a0e1a404SHong Zhang 
283a0e1a404SHong Zhang   PetscFunctionBegin;
284a0e1a404SHong Zhang   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
285595208bcSHong Zhang   ierr = MatMissingDiagonal_SeqBAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
286a0e1a404SHong Zhang 
28713f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&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 
292*f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
293*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
294f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
295ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr);
296a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
297a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
298a0e1a404SHong Zhang   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
299a0e1a404SHong Zhang 
300a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
301a0e1a404SHong Zhang   bi = b->i;
302a0e1a404SHong Zhang   bj = b->j;
303a0e1a404SHong Zhang   bv = b->a;
304a0e1a404SHong Zhang 
305a0e1a404SHong Zhang   bi[0] = 0;
306a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
307a0e1a404SHong Zhang     aj = a->j + a->diag[i];
308a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
309a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++){
310a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
311a0e1a404SHong Zhang       for (k=0; k<bs2; k++){
312a0e1a404SHong Zhang         *bv = *av; bv++; av++;
313a0e1a404SHong Zhang       }
314a0e1a404SHong Zhang     }
315a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
316a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
317a0e1a404SHong Zhang   }
318a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
319a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
321a0e1a404SHong Zhang 
322ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
3238ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
324c3d102feSKris Buschelman   } else {
325a0e1a404SHong Zhang     *newmat = B;
326c3d102feSKris Buschelman   }
327a0e1a404SHong Zhang 
328a0e1a404SHong Zhang   PetscFunctionReturn(0);
329a0e1a404SHong Zhang }
330a0e1a404SHong Zhang EXTERN_C_END
331