xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
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 EXTERN_C_BEGIN
759557b74SHong Zhang #undef __FUNCT__
84124e67cSShri Abhyankar #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqAIJ"
97087cfbeSBarry Smith PetscErrorCode  MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
104e5e7fe4SHong Zhang {
114e5e7fe4SHong Zhang   Mat            B;
124e5e7fe4SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
134e5e7fe4SHong Zhang   Mat_SeqAIJ     *b;
14dfbe8321SBarry Smith   PetscErrorCode ierr;
15d0f46423SBarry 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;
1601be0148SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0;
17dd6ea824SBarry Smith   MatScalar      *av,*bv;
184e5e7fe4SHong Zhang 
194e5e7fe4SHong Zhang   PetscFunctionBegin;
204e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
2174ed9c26SBarry Smith   ierr = PetscMalloc2(m,PetscInt,&rowlengths,m+1,PetscInt,&rowstart);CHKERRQ(ierr);
22a7a3a9ebSHong Zhang 
23a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
244e5e7fe4SHong Zhang   aj = a->j;
25a7a3a9ebSHong Zhang   k  = 0;
26a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) {
274e5e7fe4SHong Zhang     nz = ai[i+1] - ai[i];
2801be0148SBarry Smith     if (nz) {
2901be0148SBarry Smith       rowlengths[k] += nz;   /* no. of upper triangular blocks */
3001be0148SBarry Smith       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
3101be0148SBarry Smith       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
32a7a3a9ebSHong Zhang         rowlengths[(*aj)*bs]++; aj++;
33a7a3a9ebSHong Zhang       }
3401be0148SBarry Smith     }
35a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
36a7a3a9ebSHong Zhang     for (j=1; j<bs; j++) {
37a7a3a9ebSHong Zhang       rowlengths[k+j] = rowlengths[k];
38a7a3a9ebSHong Zhang     }
39a7a3a9ebSHong Zhang     k += bs;
40a7a3a9ebSHong Zhang     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
414e5e7fe4SHong Zhang   }
424e5e7fe4SHong Zhang 
437adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
44f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
456d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
46f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
474e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
48*26fbe8dcSKarl Rupp 
49d0f46423SBarry Smith   B->rmap->bs = A->rmap->bs;
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;
913bededecSBarry Smith           *(bv + rowstart[itmp]) = *(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 
111ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1128ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
113c3d102feSKris Buschelman   } else {
1144e5e7fe4SHong Zhang     *newmat = B;
115c3d102feSKris Buschelman   }
1164e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1174e5e7fe4SHong Zhang }
118be1d678aSKris Buschelman EXTERN_C_END
119be1d678aSKris Buschelman 
120be1d678aSKris Buschelman EXTERN_C_BEGIN
1214e5e7fe4SHong Zhang #undef __FUNCT__
12259557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
1230adebc6cSBarry Smith PetscErrorCode  MatConvert_SeqAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1240adebc6cSBarry Smith {
125676c34cdSKris Buschelman   Mat            B;
12659557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
127861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
128dfbe8321SBarry Smith   PetscErrorCode ierr;
129d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
130dd6ea824SBarry Smith   MatScalar      *av,*bv;
13159557b74SHong Zhang 
13259557b74SHong Zhang   PetscFunctionBegin;
133701b2c3fSHong Zhang   if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
134e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
13559557b74SHong Zhang 
13613f74950SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
13759557b74SHong Zhang   for (i=0; i<m; i++) {
13859557b74SHong Zhang     rowlengths[i] = ai[i+1] - a->diag[i];
13959557b74SHong Zhang   }
1407adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
141f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1426d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
143ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
14459557b74SHong Zhang 
1454e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);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__
1797af9e4fdSHong Zhang #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqBAIJ"
1807087cfbeSBarry Smith PetscErrorCode  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;
186d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
18774ee4d9fSHong Zhang   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
188dd6ea824SBarry Smith   MatScalar      *av,*bv;
189a0e1a404SHong Zhang 
190a0e1a404SHong Zhang   PetscFunctionBegin;
191a0e1a404SHong Zhang   /* compute browlengths of newmat */
19274ed9c26SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);CHKERRQ(ierr);
193a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
194a0e1a404SHong Zhang   aj = a->j;
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 
2047adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
205f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2066d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
207f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
2084e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
209a0e1a404SHong Zhang 
210a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
211a0e1a404SHong Zhang   bi = b->i;
212a0e1a404SHong Zhang   bj = b->j;
213a0e1a404SHong Zhang   bv = b->a;
214a0e1a404SHong Zhang 
215a0e1a404SHong Zhang   /* set b->i */
216a0e1a404SHong Zhang   bi[0] = 0;
217a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
218a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
219a0e1a404SHong Zhang     bi[i+1]      = bi[i] + browlengths[i];
220a0e1a404SHong Zhang     browstart[i] = bi[i];
221a0e1a404SHong Zhang   }
222e32f2f54SBarry 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);
223a0e1a404SHong Zhang 
224a0e1a404SHong Zhang   /* set b->j and b->a */
225a0e1a404SHong Zhang   aj = a->j; av = a->a;
226a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
227a0e1a404SHong Zhang     /* diagonal block */
228a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
229*26fbe8dcSKarl Rupp 
230a0e1a404SHong Zhang     itmp = bs2*browstart[i];
231a0e1a404SHong Zhang     for (k=0; k<bs2; k++) {
232a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
233a0e1a404SHong Zhang     }
234a0e1a404SHong Zhang     browstart[i]++;
235a0e1a404SHong Zhang 
236a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
237a0e1a404SHong Zhang     while (nz--) {
23874ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
23974ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
240*26fbe8dcSKarl Rupp 
24174ee4d9fSHong Zhang       itmp = bs2*browstart[*aj];  /* row index */
24274ee4d9fSHong Zhang       for (col=0; col<bs; col++) {
24374ee4d9fSHong Zhang         k = col;
24474ee4d9fSHong Zhang         for (row=0; row<bs; row++) {
24574ee4d9fSHong Zhang           bv[itmp + col*bs+row] = av[k]; k+=bs;
24674ee4d9fSHong Zhang         }
247a0e1a404SHong Zhang       }
248a0e1a404SHong Zhang       browstart[*aj]++;
249a0e1a404SHong Zhang 
250a0e1a404SHong Zhang       /* upper triangular blocks */
251a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
252*26fbe8dcSKarl Rupp 
253a0e1a404SHong Zhang       itmp = bs2*browstart[i];
254a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
25574ee4d9fSHong Zhang         bv[itmp + k] = av[k];
256a0e1a404SHong Zhang       }
25774ee4d9fSHong Zhang       av += bs2;
258a0e1a404SHong Zhang       browstart[i]++;
259a0e1a404SHong Zhang     }
260a0e1a404SHong Zhang   }
26174ed9c26SBarry Smith   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
262a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264a0e1a404SHong Zhang 
265ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
2668ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
267c3d102feSKris Buschelman   } else {
268a0e1a404SHong Zhang     *newmat = B;
269c3d102feSKris Buschelman   }
270a0e1a404SHong Zhang   PetscFunctionReturn(0);
271a0e1a404SHong Zhang }
272be1d678aSKris Buschelman EXTERN_C_END
273be1d678aSKris Buschelman 
274be1d678aSKris Buschelman EXTERN_C_BEGIN
275a0e1a404SHong Zhang #undef __FUNCT__
276a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
2777087cfbeSBarry Smith PetscErrorCode  MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
278a0e1a404SHong Zhang {
279a0e1a404SHong Zhang   Mat            B;
280a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
281a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
282dfbe8321SBarry Smith   PetscErrorCode ierr;
283d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
284d0f46423SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
285dd6ea824SBarry Smith   MatScalar      *av,*bv;
286ace3abfcSBarry Smith   PetscBool      flg;
287a0e1a404SHong Zhang 
288a0e1a404SHong Zhang   PetscFunctionBegin;
289701b2c3fSHong Zhang   if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
290e32f2f54SBarry Smith   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2912af78befSBarry Smith   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
292e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
293a0e1a404SHong Zhang 
29413f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
295a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
296a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
297a0e1a404SHong Zhang   }
298a0e1a404SHong Zhang 
2997adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
300f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
3016d0a4a0eSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
302ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr);
3034e0d8c25SBarry Smith   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
304a0e1a404SHong Zhang 
305a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
306a0e1a404SHong Zhang   bi = b->i;
307a0e1a404SHong Zhang   bj = b->j;
308a0e1a404SHong Zhang   bv = b->a;
309a0e1a404SHong Zhang 
310a0e1a404SHong Zhang   bi[0] = 0;
311a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
312a0e1a404SHong Zhang     aj = a->j + a->diag[i];
313a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
314a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++) {
315a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
316a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
317a0e1a404SHong Zhang         *bv = *av; bv++; av++;
318a0e1a404SHong Zhang       }
319a0e1a404SHong Zhang     }
320a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
321a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
322a0e1a404SHong Zhang   }
323a0e1a404SHong Zhang   ierr = PetscFree(browlengths);CHKERRQ(ierr);
324a0e1a404SHong Zhang   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
325a0e1a404SHong Zhang   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
326a0e1a404SHong Zhang 
327ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
3288ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
329c3d102feSKris Buschelman   } else {
330a0e1a404SHong Zhang     *newmat = B;
331c3d102feSKris Buschelman   }
332a0e1a404SHong Zhang   PetscFunctionReturn(0);
333a0e1a404SHong Zhang }
334a0e1a404SHong Zhang EXTERN_C_END
335