xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
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;
11d0f46423SBarry 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;
1201be0148SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0;
13dd6ea824SBarry Smith   MatScalar      *av,*bv;
14eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
15*b94d7dedSBarry Smith   const int      aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
16eb1ec7c1SStefano Zampini #else
17eb1ec7c1SStefano Zampini   const int      aconj = 0;
18eb1ec7c1SStefano Zampini #endif
194e5e7fe4SHong Zhang 
204e5e7fe4SHong Zhang   PetscFunctionBegin;
214e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m,&rowlengths,m+1,&rowstart));
23a7a3a9ebSHong Zhang 
24a7a3a9ebSHong Zhang   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
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;
404e5e7fe4SHong Zhang   }
414e5e7fe4SHong Zhang 
42bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
439566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
449566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,m,n,m,n));
459566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATSEQAIJ));
469566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B,0,rowlengths));
479566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(B,A->rmap->bs));
48bd019fc1SStefano Zampini   } else B = *newmat;
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   }
64aed4548fSBarry Smith   PetscCheck(bi[mbs] == 2*a->nz - diagcnt,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %" PetscInt_FMT " != 2*a->nz-diagcnt: %" PetscInt_FMT,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;
90eb1ec7c1SStefano Zampini           *(bv + rowstart[itmp]) = aconj ? PetscConj(*(av+j*bs+k)) : *(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   }
1069566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rowlengths,rowstart));
1079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1094e5e7fe4SHong Zhang 
110511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
1119566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
112c3d102feSKris Buschelman   } else {
1134e5e7fe4SHong Zhang     *newmat = B;
114c3d102feSKris Buschelman   }
1154e5e7fe4SHong Zhang   PetscFunctionReturn(0);
1164e5e7fe4SHong Zhang }
117be1d678aSKris Buschelman 
118cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1190adebc6cSBarry Smith {
120676c34cdSKris Buschelman   Mat            B;
12159557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
122861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
123ae8d29abSPierre Jolivet   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths,bs=PetscAbs(A->rmap->bs);
124dd6ea824SBarry Smith   MatScalar      *av,*bv;
125b05258aeSStefano Zampini   PetscBool      miss = PETSC_FALSE;
12659557b74SHong Zhang 
12759557b74SHong Zhang   PetscFunctionBegin;
128b05258aeSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
129*b94d7dedSBarry Smith   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
130b05258aeSStefano Zampini #else
131*b94d7dedSBarry Smith   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be either symmetric or hermitian. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE) and/or MatSetOption(mat,MAT_HERMITIAN,PETSC_TRUE)");
132b05258aeSStefano Zampini #endif
13308401ef6SPierre Jolivet   PetscCheck(n == m,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
13459557b74SHong Zhang 
1359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m/bs,&rowlengths));
136ae8d29abSPierre Jolivet   for (i=0; i<m/bs; i++) {
137b05258aeSStefano Zampini     if (a->diag[i*bs] == ai[i*bs+1]) { /* missing diagonal */
138b05258aeSStefano Zampini       rowlengths[i] = (ai[i*bs+1] - ai[i*bs])/bs; /* allocate some extra space */
139b05258aeSStefano Zampini       miss = PETSC_TRUE;
140b05258aeSStefano Zampini     } else {
141ae8d29abSPierre Jolivet       rowlengths[i] = (ai[i*bs+1] - a->diag[i*bs])/bs;
14259557b74SHong Zhang     }
143b05258aeSStefano Zampini   }
144bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
1459566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
1469566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,m,n,m,n));
1479566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATSEQSBAIJ));
1489566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B,bs,0,rowlengths));
149bd019fc1SStefano Zampini   } else B = *newmat;
15059557b74SHong Zhang 
151b05258aeSStefano Zampini   if (bs == 1 && !miss) {
152676c34cdSKris Buschelman     b  = (Mat_SeqSBAIJ*)(B->data);
153861ba921SHong Zhang     bi = b->i;
154861ba921SHong Zhang     bj = b->j;
155861ba921SHong Zhang     bv = b->a;
156861ba921SHong Zhang 
157861ba921SHong Zhang     bi[0] = 0;
15859557b74SHong Zhang     for (i=0; i<m; i++) {
15959557b74SHong Zhang       aj = a->j + a->diag[i];
16059557b74SHong Zhang       av = a->a + a->diag[i];
161861ba921SHong Zhang       for (j=0; j<rowlengths[i]; j++) {
162861ba921SHong Zhang         *bj = *aj; bj++; aj++;
163861ba921SHong Zhang         *bv = *av; bv++; av++;
164861ba921SHong Zhang       }
165861ba921SHong Zhang       bi[i+1]    = bi[i] + rowlengths[i];
166861ba921SHong Zhang       b->ilen[i] = rowlengths[i];
16759557b74SHong Zhang     }
1689566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
170ae8d29abSPierre Jolivet   } else {
1719566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
172ae8d29abSPierre 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 */
173ae8d29abSPierre Jolivet     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
174ae8d29abSPierre Jolivet     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
1759566063dSJacob Faibussowitsch     PetscCall(MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&B));
176ae8d29abSPierre Jolivet   }
1779566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowlengths));
178511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
1799566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
180ae8d29abSPierre Jolivet   } else *newmat = B;
18159557b74SHong Zhang   PetscFunctionReturn(0);
18259557b74SHong Zhang }
18359557b74SHong Zhang 
184cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
185a0e1a404SHong Zhang {
186a0e1a404SHong Zhang   Mat            B;
187a0e1a404SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
188a0e1a404SHong Zhang   Mat_SeqBAIJ    *b;
189d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
19074ee4d9fSHong Zhang   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
191dd6ea824SBarry Smith   MatScalar      *av,*bv;
192eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
193*b94d7dedSBarry Smith   const int      aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
194eb1ec7c1SStefano Zampini #else
195eb1ec7c1SStefano Zampini   const int      aconj = 0;
196eb1ec7c1SStefano Zampini #endif
197a0e1a404SHong Zhang 
198a0e1a404SHong Zhang   PetscFunctionBegin;
199a0e1a404SHong Zhang   /* compute browlengths of newmat */
2009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs,&browlengths,mbs,&browstart));
201a0e1a404SHong Zhang   for (i=0; i<mbs; i++) browlengths[i] = 0;
202a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
203a0e1a404SHong Zhang     nz = ai[i+1] - ai[i];
204a0e1a404SHong Zhang     aj++; /* skip diagonal */
205a0e1a404SHong Zhang     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
206a0e1a404SHong Zhang       browlengths[*aj]++; aj++;
207a0e1a404SHong Zhang     }
208a0e1a404SHong Zhang     browlengths[i] += nz;   /* no. of upper triangular blocks */
209a0e1a404SHong Zhang   }
210a0e1a404SHong Zhang 
211bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2129566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
2139566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,m,n,m,n));
2149566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATSEQBAIJ));
2159566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetPreallocation(B,bs,0,browlengths));
216bd019fc1SStefano Zampini   } else B = *newmat;
217a0e1a404SHong Zhang 
218a0e1a404SHong Zhang   b  = (Mat_SeqBAIJ*)(B->data);
219a0e1a404SHong Zhang   bi = b->i;
220a0e1a404SHong Zhang   bj = b->j;
221a0e1a404SHong Zhang   bv = b->a;
222a0e1a404SHong Zhang 
223a0e1a404SHong Zhang   /* set b->i */
224a0e1a404SHong Zhang   bi[0] = 0;
225a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
226a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
227a0e1a404SHong Zhang     bi[i+1]      = bi[i] + browlengths[i];
228a0e1a404SHong Zhang     browstart[i] = bi[i];
229a0e1a404SHong Zhang   }
230aed4548fSBarry Smith   PetscCheck(bi[mbs] == 2*a->nz - mbs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %" PetscInt_FMT " != 2*a->nz - mbs: %" PetscInt_FMT,bi[mbs],2*a->nz - mbs);
231a0e1a404SHong Zhang 
232a0e1a404SHong Zhang   /* set b->j and b->a */
233a0e1a404SHong Zhang   aj = a->j; av = a->a;
234a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
235a0e1a404SHong Zhang     /* diagonal block */
236a0e1a404SHong Zhang     *(bj + browstart[i]) = *aj; aj++;
23726fbe8dcSKarl Rupp 
238a0e1a404SHong Zhang     itmp = bs2*browstart[i];
239a0e1a404SHong Zhang     for (k=0; k<bs2; k++) {
240a0e1a404SHong Zhang       *(bv + itmp + k) = *av; av++;
241a0e1a404SHong Zhang     }
242a0e1a404SHong Zhang     browstart[i]++;
243a0e1a404SHong Zhang 
244a0e1a404SHong Zhang     nz = ai[i+1] - ai[i] -1;
245a0e1a404SHong Zhang     while (nz--) {
24674ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
24774ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
24826fbe8dcSKarl Rupp 
24974ee4d9fSHong Zhang       itmp = bs2*browstart[*aj];  /* row index */
25074ee4d9fSHong Zhang       for (col=0; col<bs; col++) {
25174ee4d9fSHong Zhang         k = col;
25274ee4d9fSHong Zhang         for (row=0; row<bs; row++) {
253eb1ec7c1SStefano Zampini           bv[itmp + col*bs+row] = aconj ? PetscConj(av[k]) : av[k];
254eb1ec7c1SStefano Zampini           k+=bs;
25574ee4d9fSHong Zhang         }
256a0e1a404SHong Zhang       }
257a0e1a404SHong Zhang       browstart[*aj]++;
258a0e1a404SHong Zhang 
259a0e1a404SHong Zhang       /* upper triangular blocks */
260a0e1a404SHong Zhang       *(bj + browstart[i]) = *aj; aj++;
26126fbe8dcSKarl Rupp 
262a0e1a404SHong Zhang       itmp = bs2*browstart[i];
263a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
26474ee4d9fSHong Zhang         bv[itmp + k] = av[k];
265a0e1a404SHong Zhang       }
26674ee4d9fSHong Zhang       av += bs2;
267a0e1a404SHong Zhang       browstart[i]++;
268a0e1a404SHong Zhang     }
269a0e1a404SHong Zhang   }
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree2(browlengths,browstart));
2719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
273a0e1a404SHong Zhang 
274511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
2759566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
276ae8d29abSPierre Jolivet   } else *newmat = B;
277a0e1a404SHong Zhang   PetscFunctionReturn(0);
278a0e1a404SHong Zhang }
279be1d678aSKris Buschelman 
280cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
281a0e1a404SHong Zhang {
282a0e1a404SHong Zhang   Mat            B;
283a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
284a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
285d0f46423SBarry Smith   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
286d0f46423SBarry Smith   PetscInt       bs =A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
287dd6ea824SBarry Smith   MatScalar      *av,*bv;
288ace3abfcSBarry Smith   PetscBool      flg;
289a0e1a404SHong Zhang 
290a0e1a404SHong Zhang   PetscFunctionBegin;
29128b400f6SJacob Faibussowitsch   PetscCheck(A->symmetric,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
29208401ef6SPierre Jolivet   PetscCheck(n == m,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
2939566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal_SeqBAIJ(A,&flg,&dd)); /* check for missing diagonals, then mark diag */
29428b400f6SJacob Faibussowitsch   PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %" PetscInt_FMT,dd);
295a0e1a404SHong Zhang 
2969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs,&browlengths));
297a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
298a0e1a404SHong Zhang     browlengths[i] = ai[i+1] - a->diag[i];
299a0e1a404SHong Zhang   }
300a0e1a404SHong Zhang 
301bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
3029566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
3039566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,m,n,m,n));
3049566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATSEQSBAIJ));
3059566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B,bs,0,browlengths));
306bd019fc1SStefano Zampini   } else B = *newmat;
307a0e1a404SHong Zhang 
308a0e1a404SHong Zhang   b  = (Mat_SeqSBAIJ*)(B->data);
309a0e1a404SHong Zhang   bi = b->i;
310a0e1a404SHong Zhang   bj = b->j;
311a0e1a404SHong Zhang   bv = b->a;
312a0e1a404SHong Zhang 
313a0e1a404SHong Zhang   bi[0] = 0;
314a0e1a404SHong Zhang   for (i=0; i<mbs; i++) {
315a0e1a404SHong Zhang     aj = a->j + a->diag[i];
316a0e1a404SHong Zhang     av = a->a + (a->diag[i])*bs2;
317a0e1a404SHong Zhang     for (j=0; j<browlengths[i]; j++) {
318a0e1a404SHong Zhang       *bj = *aj; bj++; aj++;
319a0e1a404SHong Zhang       for (k=0; k<bs2; k++) {
320a0e1a404SHong Zhang         *bv = *av; bv++; av++;
321a0e1a404SHong Zhang       }
322a0e1a404SHong Zhang     }
323a0e1a404SHong Zhang     bi[i+1]    = bi[i] + browlengths[i];
324a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
325a0e1a404SHong Zhang   }
3269566063dSJacob Faibussowitsch   PetscCall(PetscFree(browlengths));
3279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
3289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
329a0e1a404SHong Zhang 
330511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
3319566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
332ae8d29abSPierre Jolivet   } else *newmat = B;
333a0e1a404SHong Zhang   PetscFunctionReturn(0);
334a0e1a404SHong Zhang }
335