xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision ae1ee55146a7ad071171b861759b85940c7e5c67)
1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
459557b74SHong Zhang 
MatConvert_SeqSBAIJ_SeqAIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)5d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
6d71ae5a4SJacob Faibussowitsch {
74e5e7fe4SHong Zhang   Mat           B;
84e5e7fe4SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
94e5e7fe4SHong Zhang   Mat_SeqAIJ   *b;
10d0f46423SBarry 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;
1101be0148SBarry Smith   PetscInt      bs = A->rmap->bs, bs2 = bs * bs, mbs = A->rmap->N / bs, diagcnt = 0;
12dd6ea824SBarry Smith   MatScalar    *av, *bv;
13*fc2fb351SPierre Jolivet   const int     aconj = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
144e5e7fe4SHong Zhang 
154e5e7fe4SHong Zhang   PetscFunctionBegin;
164e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m, &rowlengths, m + 1, &rowstart));
18a7a3a9ebSHong Zhang 
19a7a3a9ebSHong Zhang   for (i = 0; i < mbs; i++) rowlengths[i * bs] = 0;
20a7a3a9ebSHong Zhang   k = 0;
21a7a3a9ebSHong Zhang   for (i = 0; i < mbs; i++) {
224e5e7fe4SHong Zhang     nz = ai[i + 1] - ai[i];
2301be0148SBarry Smith     if (nz) {
2401be0148SBarry Smith       rowlengths[k] += nz; /* no. of upper triangular blocks */
259371c9d4SSatish Balay       if (*aj == i) {
269371c9d4SSatish Balay         aj++;
279371c9d4SSatish Balay         diagcnt++;
289371c9d4SSatish Balay         nz--;
299371c9d4SSatish Balay       } /* skip diagonal */
3001be0148SBarry Smith       for (j = 0; j < nz; j++) { /* no. of lower triangular blocks */
319371c9d4SSatish Balay         rowlengths[(*aj) * bs]++;
329371c9d4SSatish Balay         aj++;
33a7a3a9ebSHong Zhang       }
3401be0148SBarry Smith     }
35a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
36ad540459SPierre Jolivet     for (j = 1; j < bs; j++) rowlengths[k + j] = rowlengths[k];
37a7a3a9ebSHong Zhang     k += bs;
384e5e7fe4SHong Zhang   }
394e5e7fe4SHong Zhang 
40bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
419566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
429566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
439566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
449566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, rowlengths));
459566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(B, A->rmap->bs));
46bd019fc1SStefano Zampini   } else B = *newmat;
474e5e7fe4SHong Zhang 
48f4f49eeaSPierre Jolivet   b  = (Mat_SeqAIJ *)B->data;
494e5e7fe4SHong Zhang   bi = b->i;
504e5e7fe4SHong Zhang   bj = b->j;
514e5e7fe4SHong Zhang   bv = b->a;
524e5e7fe4SHong Zhang 
534e5e7fe4SHong Zhang   /* set b->i */
549371c9d4SSatish Balay   bi[0]       = 0;
559371c9d4SSatish Balay   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   }
63aed4548fSBarry 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);
644e5e7fe4SHong Zhang 
654e5e7fe4SHong Zhang   /* set b->j and b->a */
669371c9d4SSatish Balay   aj = a->j;
679371c9d4SSatish Balay   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       }
819371c9d4SSatish Balay       aj++;
829371c9d4SSatish Balay       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;
91eb1ec7c1SStefano 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       }
1049371c9d4SSatish Balay       aj++;
1059371c9d4SSatish Balay       av += bs2;
1064e5e7fe4SHong Zhang     }
1074e5e7fe4SHong Zhang   }
1089566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rowlengths, rowstart));
1099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1109566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1114e5e7fe4SHong Zhang 
112511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
1139566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
114c3d102feSKris Buschelman   } else {
1154e5e7fe4SHong Zhang     *newmat = B;
116c3d102feSKris Buschelman   }
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1184e5e7fe4SHong Zhang }
119be1d678aSKris Buschelman 
1205a2b941aSBarry Smith #include <../src/mat/impls/aij/seq/aij.h>
1215a2b941aSBarry Smith 
MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(Mat A,PetscInt ** nnz)1225a2b941aSBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(Mat A, PetscInt **nnz)
1235a2b941aSBarry Smith {
1245a2b941aSBarry Smith   Mat_SeqAIJ     *Aa = (Mat_SeqAIJ *)A->data;
12558b7e2c1SStefano Zampini   PetscInt        m, n, bs = A->rmap->bs;
126f85a0629SBarry Smith   const PetscInt *ai = Aa->i, *aj = Aa->j;
1275a2b941aSBarry Smith 
1285a2b941aSBarry Smith   PetscFunctionBegin;
1295a2b941aSBarry Smith   PetscCall(MatGetSize(A, &m, &n));
1305a2b941aSBarry Smith 
131f85a0629SBarry Smith   if (bs == 1) {
132421480d9SBarry Smith     const PetscInt *adiag;
133f85a0629SBarry Smith 
134421480d9SBarry Smith     PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
135f85a0629SBarry Smith     PetscCall(PetscMalloc1(m, nnz));
136f85a0629SBarry Smith     for (PetscInt i = 0; i < m; i++) {
137f85a0629SBarry Smith       if (adiag[i] == ai[i + 1]) {
138f85a0629SBarry Smith         (*nnz)[i] = 0;
139f85a0629SBarry Smith         for (PetscInt j = ai[i]; j < ai[i + 1]; j++) (*nnz)[i] += (aj[j] > i);
140f85a0629SBarry Smith       } else (*nnz)[i] = ai[i + 1] - adiag[i];
141f85a0629SBarry Smith     }
142f85a0629SBarry Smith   } else {
143f85a0629SBarry Smith     PetscHSetIJ    ht;
144f85a0629SBarry Smith     PetscHashIJKey key;
145f85a0629SBarry Smith     PetscBool      missing;
146f85a0629SBarry Smith 
147f85a0629SBarry Smith     PetscCall(PetscHSetIJCreate(&ht));
148f85a0629SBarry Smith     PetscCall(PetscCalloc1(m / bs, nnz));
149f85a0629SBarry Smith     for (PetscInt i = 0; i < m; i++) {
150f85a0629SBarry Smith       key.i = i / bs;
151f85a0629SBarry Smith       for (PetscInt k = ai[i]; k < ai[i + 1]; k++) {
152f85a0629SBarry Smith         key.j = aj[k] / bs;
153f85a0629SBarry Smith         if (key.j < key.i) continue;
154f85a0629SBarry Smith         PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
155f85a0629SBarry Smith         if (missing) (*nnz)[key.i]++;
1565a2b941aSBarry Smith       }
1575a2b941aSBarry Smith     }
158f85a0629SBarry Smith     PetscCall(PetscHSetIJDestroy(&ht));
1595a2b941aSBarry Smith   }
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1615a2b941aSBarry Smith }
1625a2b941aSBarry Smith 
MatConvert_SeqAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)163d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
164d71ae5a4SJacob Faibussowitsch {
165676c34cdSKris Buschelman   Mat             B;
16659557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
167861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
16858b7e2c1SStefano Zampini   PetscInt       *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->N, i, j, *bi, *bj, *rowlengths, bs = A->rmap->bs;
169dd6ea824SBarry Smith   MatScalar      *av, *bv;
170b05258aeSStefano Zampini   PetscBool       miss = PETSC_FALSE;
171421480d9SBarry Smith   const PetscInt *adiag;
17259557b74SHong Zhang 
17359557b74SHong Zhang   PetscFunctionBegin;
174b05258aeSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
175b94d7dedSBarry Smith   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
176b05258aeSStefano Zampini #else
177b0c98d1dSPierre Jolivet   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)");
178b05258aeSStefano Zampini #endif
17908401ef6SPierre Jolivet   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square");
180421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
1815a2b941aSBarry Smith   if (bs == 1) {
1825a2b941aSBarry Smith     PetscCall(PetscMalloc1(m, &rowlengths));
1835a2b941aSBarry Smith     for (i = 0; i < m; i++) {
184421480d9SBarry Smith       if (adiag[i] == ai[i + 1]) {               /* missing diagonal */
1855a2b941aSBarry Smith         rowlengths[i] = (ai[i + 1] - ai[i]) + 1; /* allocate some extra space */
186b05258aeSStefano Zampini         miss          = PETSC_TRUE;
187b05258aeSStefano Zampini       } else {
188421480d9SBarry Smith         rowlengths[i] = (ai[i + 1] - adiag[i]);
18959557b74SHong Zhang       }
190b05258aeSStefano Zampini     }
1915a2b941aSBarry Smith   } else PetscCall(MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(A, &rowlengths));
192bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
1939566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1949566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
1959566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSBAIJ));
1969566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, rowlengths));
197bd019fc1SStefano Zampini   } else B = *newmat;
19859557b74SHong Zhang 
199b05258aeSStefano Zampini   if (bs == 1 && !miss) {
200f4f49eeaSPierre Jolivet     b  = (Mat_SeqSBAIJ *)B->data;
201861ba921SHong Zhang     bi = b->i;
202861ba921SHong Zhang     bj = b->j;
203861ba921SHong Zhang     bv = b->a;
204861ba921SHong Zhang 
205861ba921SHong Zhang     bi[0] = 0;
20659557b74SHong Zhang     for (i = 0; i < m; i++) {
207421480d9SBarry Smith       aj = a->j + adiag[i];
208421480d9SBarry Smith       av = a->a + adiag[i];
209861ba921SHong Zhang       for (j = 0; j < rowlengths[i]; j++) {
2109371c9d4SSatish Balay         *bj = *aj;
2119371c9d4SSatish Balay         bj++;
2129371c9d4SSatish Balay         aj++;
2139371c9d4SSatish Balay         *bv = *av;
2149371c9d4SSatish Balay         bv++;
2159371c9d4SSatish Balay         av++;
216861ba921SHong Zhang       }
217861ba921SHong Zhang       bi[i + 1]  = bi[i] + rowlengths[i];
218861ba921SHong Zhang       b->ilen[i] = rowlengths[i];
21959557b74SHong Zhang     }
2209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
222ae8d29abSPierre Jolivet   } else {
2239566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
224ae8d29abSPierre 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 */
225ae8d29abSPierre Jolivet     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
226ae8d29abSPierre Jolivet     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
2279566063dSJacob Faibussowitsch     PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B));
228ae8d29abSPierre Jolivet   }
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowlengths));
230ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
231ac530a7eSPierre Jolivet   else *newmat = B;
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23359557b74SHong Zhang }
23459557b74SHong Zhang 
MatConvert_SeqSBAIJ_SeqBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)235d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
236d71ae5a4SJacob Faibussowitsch {
237a0e1a404SHong Zhang   Mat           B;
238a0e1a404SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
239a0e1a404SHong Zhang   Mat_SeqBAIJ  *b;
240421480d9SBarry Smith   PetscInt     *ai = a->i, *aj = a->j, m = A->rmap->N, n = A->cmap->n, *bi, *bj, *browlengths, nz, *browstart, itmp;
241421480d9SBarry Smith   PetscInt      bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs;
242dd6ea824SBarry Smith   MatScalar    *av, *bv;
243*fc2fb351SPierre Jolivet   const int     aconj = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
244a0e1a404SHong Zhang 
245a0e1a404SHong Zhang   PetscFunctionBegin;
246a0e1a404SHong Zhang   /* compute browlengths of newmat */
2479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &browlengths, mbs, &browstart));
248421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) browlengths[i] = 0;
249421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
250a0e1a404SHong Zhang     nz = ai[i + 1] - ai[i];
251a0e1a404SHong Zhang     aj++;                               /* skip diagonal */
252421480d9SBarry Smith     for (PetscInt k = 1; k < nz; k++) { /* no. of lower triangular blocks */
2539371c9d4SSatish Balay       browlengths[*aj]++;
2549371c9d4SSatish Balay       aj++;
255a0e1a404SHong Zhang     }
256a0e1a404SHong Zhang     browlengths[i] += nz; /* no. of upper triangular blocks */
257a0e1a404SHong Zhang   }
258a0e1a404SHong Zhang 
259bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2609566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2619566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
2629566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQBAIJ));
2639566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetPreallocation(B, bs, 0, browlengths));
264bd019fc1SStefano Zampini   } else B = *newmat;
265a0e1a404SHong Zhang 
266f4f49eeaSPierre Jolivet   b  = (Mat_SeqBAIJ *)B->data;
267a0e1a404SHong Zhang   bi = b->i;
268a0e1a404SHong Zhang   bj = b->j;
269a0e1a404SHong Zhang   bv = b->a;
270a0e1a404SHong Zhang 
271a0e1a404SHong Zhang   /* set b->i */
272a0e1a404SHong Zhang   bi[0] = 0;
273421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
274a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
275a0e1a404SHong Zhang     bi[i + 1]    = bi[i] + browlengths[i];
276a0e1a404SHong Zhang     browstart[i] = bi[i];
277a0e1a404SHong Zhang   }
278aed4548fSBarry 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);
279a0e1a404SHong Zhang 
280a0e1a404SHong Zhang   /* set b->j and b->a */
2819371c9d4SSatish Balay   aj = a->j;
2829371c9d4SSatish Balay   av = a->a;
283421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
284a0e1a404SHong Zhang     /* diagonal block */
2859371c9d4SSatish Balay     *(bj + browstart[i]) = *aj;
2869371c9d4SSatish Balay     aj++;
28726fbe8dcSKarl Rupp 
288a0e1a404SHong Zhang     itmp = bs2 * browstart[i];
289421480d9SBarry Smith     for (PetscInt k = 0; k < bs2; k++) {
2909371c9d4SSatish Balay       *(bv + itmp + k) = *av;
2919371c9d4SSatish Balay       av++;
292a0e1a404SHong Zhang     }
293a0e1a404SHong Zhang     browstart[i]++;
294a0e1a404SHong Zhang 
295a0e1a404SHong Zhang     nz = ai[i + 1] - ai[i] - 1;
296a0e1a404SHong Zhang     while (nz--) {
29774ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
29874ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
29926fbe8dcSKarl Rupp 
30074ee4d9fSHong Zhang       itmp = bs2 * browstart[*aj]; /* row index */
301421480d9SBarry Smith       for (PetscInt col = 0; col < bs; col++) {
302421480d9SBarry Smith         PetscInt k = col;
303421480d9SBarry Smith         for (PetscInt row = 0; row < bs; row++) {
304eb1ec7c1SStefano Zampini           bv[itmp + col * bs + row] = aconj ? PetscConj(av[k]) : av[k];
305eb1ec7c1SStefano Zampini           k += bs;
30674ee4d9fSHong Zhang         }
307a0e1a404SHong Zhang       }
308a0e1a404SHong Zhang       browstart[*aj]++;
309a0e1a404SHong Zhang 
310a0e1a404SHong Zhang       /* upper triangular blocks */
3119371c9d4SSatish Balay       *(bj + browstart[i]) = *aj;
3129371c9d4SSatish Balay       aj++;
31326fbe8dcSKarl Rupp 
314a0e1a404SHong Zhang       itmp = bs2 * browstart[i];
315421480d9SBarry Smith       for (PetscInt k = 0; k < bs2; k++) bv[itmp + k] = av[k];
31674ee4d9fSHong Zhang       av += bs2;
317a0e1a404SHong Zhang       browstart[i]++;
318a0e1a404SHong Zhang     }
319a0e1a404SHong Zhang   }
3209566063dSJacob Faibussowitsch   PetscCall(PetscFree2(browlengths, browstart));
3219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
323a0e1a404SHong Zhang 
324ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
325ac530a7eSPierre Jolivet   else *newmat = B;
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327a0e1a404SHong Zhang }
328be1d678aSKris Buschelman 
MatConvert_SeqBAIJ_SeqSBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)329d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
330d71ae5a4SJacob Faibussowitsch {
331a0e1a404SHong Zhang   Mat             B;
332a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data;
333a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
334421480d9SBarry Smith   PetscInt       *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->n, k, *bi, *bj, *browlengths;
335421480d9SBarry Smith   PetscInt        bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs;
336dd6ea824SBarry Smith   MatScalar      *av, *bv;
337421480d9SBarry Smith   const PetscInt *adiag;
338a0e1a404SHong Zhang 
339a0e1a404SHong Zhang   PetscFunctionBegin;
340b0c98d1dSPierre Jolivet   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
34108401ef6SPierre Jolivet   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square");
342421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &adiag, NULL));
3439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs, &browlengths));
344421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) browlengths[i] = ai[i + 1] - adiag[i];
345a0e1a404SHong Zhang 
346bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
3479566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
3489566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
3499566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSBAIJ));
3509566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, browlengths));
351bd019fc1SStefano Zampini   } else B = *newmat;
352a0e1a404SHong Zhang 
353f4f49eeaSPierre Jolivet   b  = (Mat_SeqSBAIJ *)B->data;
354a0e1a404SHong Zhang   bi = b->i;
355a0e1a404SHong Zhang   bj = b->j;
356a0e1a404SHong Zhang   bv = b->a;
357a0e1a404SHong Zhang 
358a0e1a404SHong Zhang   bi[0] = 0;
359421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
360421480d9SBarry Smith     aj = a->j + adiag[i];
361421480d9SBarry Smith     av = a->a + (adiag[i]) * bs2;
362421480d9SBarry Smith     for (PetscInt j = 0; j < browlengths[i]; j++) {
3639371c9d4SSatish Balay       *bj = *aj;
3649371c9d4SSatish Balay       bj++;
3659371c9d4SSatish Balay       aj++;
366a0e1a404SHong Zhang       for (k = 0; k < bs2; k++) {
3679371c9d4SSatish Balay         *bv = *av;
3689371c9d4SSatish Balay         bv++;
3699371c9d4SSatish Balay         av++;
370a0e1a404SHong Zhang       }
371a0e1a404SHong Zhang     }
372a0e1a404SHong Zhang     bi[i + 1]  = bi[i] + browlengths[i];
373a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
374a0e1a404SHong Zhang   }
3759566063dSJacob Faibussowitsch   PetscCall(PetscFree(browlengths));
3769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
378a0e1a404SHong Zhang 
379ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
380ac530a7eSPierre Jolivet   else *newmat = B;
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
382a0e1a404SHong Zhang }
383