xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision b0c98d1d8bc8fbb369fd6b04fbfd2a9276aa7d86)
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 
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;
13eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
14b94d7dedSBarry Smith   const int aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
15eb1ec7c1SStefano Zampini #else
16eb1ec7c1SStefano Zampini   const int aconj = 0;
17eb1ec7c1SStefano Zampini #endif
184e5e7fe4SHong Zhang 
194e5e7fe4SHong Zhang   PetscFunctionBegin;
204e5e7fe4SHong Zhang   /* compute rowlengths of newmat */
219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m, &rowlengths, m + 1, &rowstart));
22a7a3a9ebSHong Zhang 
23a7a3a9ebSHong Zhang   for (i = 0; i < mbs; i++) rowlengths[i * bs] = 0;
24a7a3a9ebSHong Zhang   k = 0;
25a7a3a9ebSHong Zhang   for (i = 0; i < mbs; i++) {
264e5e7fe4SHong Zhang     nz = ai[i + 1] - ai[i];
2701be0148SBarry Smith     if (nz) {
2801be0148SBarry Smith       rowlengths[k] += nz; /* no. of upper triangular blocks */
299371c9d4SSatish Balay       if (*aj == i) {
309371c9d4SSatish Balay         aj++;
319371c9d4SSatish Balay         diagcnt++;
329371c9d4SSatish Balay         nz--;
339371c9d4SSatish Balay       } /* skip diagonal */
3401be0148SBarry Smith       for (j = 0; j < nz; j++) { /* no. of lower triangular blocks */
359371c9d4SSatish Balay         rowlengths[(*aj) * bs]++;
369371c9d4SSatish Balay         aj++;
37a7a3a9ebSHong Zhang       }
3801be0148SBarry Smith     }
39a7a3a9ebSHong Zhang     rowlengths[k] *= bs;
40ad540459SPierre Jolivet     for (j = 1; j < bs; j++) rowlengths[k + j] = rowlengths[k];
41a7a3a9ebSHong Zhang     k += bs;
424e5e7fe4SHong Zhang   }
434e5e7fe4SHong Zhang 
44bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
459566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
469566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
479566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
489566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, rowlengths));
499566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(B, A->rmap->bs));
50bd019fc1SStefano Zampini   } else B = *newmat;
514e5e7fe4SHong Zhang 
52f4f49eeaSPierre Jolivet   b  = (Mat_SeqAIJ *)B->data;
534e5e7fe4SHong Zhang   bi = b->i;
544e5e7fe4SHong Zhang   bj = b->j;
554e5e7fe4SHong Zhang   bv = b->a;
564e5e7fe4SHong Zhang 
574e5e7fe4SHong Zhang   /* set b->i */
589371c9d4SSatish Balay   bi[0]       = 0;
599371c9d4SSatish Balay   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   }
67aed4548fSBarry 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);
684e5e7fe4SHong Zhang 
694e5e7fe4SHong Zhang   /* set b->j and b->a */
709371c9d4SSatish Balay   aj = a->j;
719371c9d4SSatish Balay   av = a->a;
72a7a3a9ebSHong Zhang   for (i = 0; i < mbs; i++) {
7301be0148SBarry Smith     nz = ai[i + 1] - ai[i];
74a7a3a9ebSHong Zhang     /* diagonal block */
7501be0148SBarry Smith     if (nz && *aj == i) {
7601be0148SBarry Smith       nz--;
77a7a3a9ebSHong Zhang       for (j = 0; j < bs; j++) { /* row i*bs+j */
78a7a3a9ebSHong Zhang         itmp = i * bs + j;
79a7a3a9ebSHong Zhang         for (k = 0; k < bs; k++) { /* col i*bs+k */
80a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj) * bs + k;
81a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av + k * bs + j);
82a7a3a9ebSHong Zhang           rowstart[itmp]++;
83a7a3a9ebSHong Zhang         }
84a7a3a9ebSHong Zhang       }
859371c9d4SSatish Balay       aj++;
869371c9d4SSatish Balay       av += bs2;
8701be0148SBarry Smith     }
88a7a3a9ebSHong Zhang 
894e5e7fe4SHong Zhang     while (nz--) {
90a7a3a9ebSHong Zhang       /* lower triangular blocks */
91a7a3a9ebSHong Zhang       for (j = 0; j < bs; j++) { /* row (*aj)*bs+j */
92a7a3a9ebSHong Zhang         itmp = (*aj) * bs + j;
93a7a3a9ebSHong Zhang         for (k = 0; k < bs; k++) { /* col i*bs+k */
94a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = i * bs + k;
95eb1ec7c1SStefano Zampini           *(bv + rowstart[itmp]) = aconj ? PetscConj(*(av + j * bs + k)) : *(av + j * bs + k);
96a7a3a9ebSHong Zhang           rowstart[itmp]++;
97a7a3a9ebSHong Zhang         }
98a7a3a9ebSHong Zhang       }
99a7a3a9ebSHong Zhang       /* upper triangular blocks */
100a7a3a9ebSHong Zhang       for (j = 0; j < bs; j++) { /* row i*bs+j */
101a7a3a9ebSHong Zhang         itmp = i * bs + j;
102a7a3a9ebSHong Zhang         for (k = 0; k < bs; k++) { /* col (*aj)*bs+k */
103a7a3a9ebSHong Zhang           *(bj + rowstart[itmp]) = (*aj) * bs + k;
104a7a3a9ebSHong Zhang           *(bv + rowstart[itmp]) = *(av + k * bs + j);
105a7a3a9ebSHong Zhang           rowstart[itmp]++;
106a7a3a9ebSHong Zhang         }
107a7a3a9ebSHong Zhang       }
1089371c9d4SSatish Balay       aj++;
1099371c9d4SSatish Balay       av += bs2;
1104e5e7fe4SHong Zhang     }
1114e5e7fe4SHong Zhang   }
1129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rowlengths, rowstart));
1139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1154e5e7fe4SHong Zhang 
116511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
1179566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
118c3d102feSKris Buschelman   } else {
1194e5e7fe4SHong Zhang     *newmat = B;
120c3d102feSKris Buschelman   }
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1224e5e7fe4SHong Zhang }
123be1d678aSKris Buschelman 
1245a2b941aSBarry Smith #include <../src/mat/impls/aij/seq/aij.h>
1255a2b941aSBarry Smith 
1265a2b941aSBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(Mat A, PetscInt **nnz)
1275a2b941aSBarry Smith {
1285a2b941aSBarry Smith   Mat_SeqAIJ     *Aa = (Mat_SeqAIJ *)A->data;
12958b7e2c1SStefano Zampini   PetscInt        m, n, bs = A->rmap->bs;
130f85a0629SBarry Smith   const PetscInt *ai = Aa->i, *aj = Aa->j;
1315a2b941aSBarry Smith 
1325a2b941aSBarry Smith   PetscFunctionBegin;
1335a2b941aSBarry Smith   PetscCall(MatGetSize(A, &m, &n));
1345a2b941aSBarry Smith 
135f85a0629SBarry Smith   if (bs == 1) {
136421480d9SBarry Smith     const PetscInt *adiag;
137f85a0629SBarry Smith 
138421480d9SBarry Smith     PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
139f85a0629SBarry Smith     PetscCall(PetscMalloc1(m, nnz));
140f85a0629SBarry Smith     for (PetscInt i = 0; i < m; i++) {
141f85a0629SBarry Smith       if (adiag[i] == ai[i + 1]) {
142f85a0629SBarry Smith         (*nnz)[i] = 0;
143f85a0629SBarry Smith         for (PetscInt j = ai[i]; j < ai[i + 1]; j++) (*nnz)[i] += (aj[j] > i);
144f85a0629SBarry Smith       } else (*nnz)[i] = ai[i + 1] - adiag[i];
145f85a0629SBarry Smith     }
146f85a0629SBarry Smith   } else {
147f85a0629SBarry Smith     PetscHSetIJ    ht;
148f85a0629SBarry Smith     PetscHashIJKey key;
149f85a0629SBarry Smith     PetscBool      missing;
150f85a0629SBarry Smith 
151f85a0629SBarry Smith     PetscCall(PetscHSetIJCreate(&ht));
152f85a0629SBarry Smith     PetscCall(PetscCalloc1(m / bs, nnz));
153f85a0629SBarry Smith     for (PetscInt i = 0; i < m; i++) {
154f85a0629SBarry Smith       key.i = i / bs;
155f85a0629SBarry Smith       for (PetscInt k = ai[i]; k < ai[i + 1]; k++) {
156f85a0629SBarry Smith         key.j = aj[k] / bs;
157f85a0629SBarry Smith         if (key.j < key.i) continue;
158f85a0629SBarry Smith         PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
159f85a0629SBarry Smith         if (missing) (*nnz)[key.i]++;
1605a2b941aSBarry Smith       }
1615a2b941aSBarry Smith     }
162f85a0629SBarry Smith     PetscCall(PetscHSetIJDestroy(&ht));
1635a2b941aSBarry Smith   }
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1655a2b941aSBarry Smith }
1665a2b941aSBarry Smith 
167d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
168d71ae5a4SJacob Faibussowitsch {
169676c34cdSKris Buschelman   Mat             B;
17059557b74SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
171861ba921SHong Zhang   Mat_SeqSBAIJ   *b;
17258b7e2c1SStefano Zampini   PetscInt       *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->N, i, j, *bi, *bj, *rowlengths, bs = A->rmap->bs;
173dd6ea824SBarry Smith   MatScalar      *av, *bv;
174b05258aeSStefano Zampini   PetscBool       miss = PETSC_FALSE;
175421480d9SBarry Smith   const PetscInt *adiag;
17659557b74SHong Zhang 
17759557b74SHong Zhang   PetscFunctionBegin;
178b05258aeSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
179b94d7dedSBarry Smith   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
180b05258aeSStefano Zampini #else
181*b0c98d1dSPierre 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)");
182b05258aeSStefano Zampini #endif
18308401ef6SPierre Jolivet   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square");
184421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
1855a2b941aSBarry Smith   if (bs == 1) {
1865a2b941aSBarry Smith     PetscCall(PetscMalloc1(m, &rowlengths));
1875a2b941aSBarry Smith     for (i = 0; i < m; i++) {
188421480d9SBarry Smith       if (adiag[i] == ai[i + 1]) {               /* missing diagonal */
1895a2b941aSBarry Smith         rowlengths[i] = (ai[i + 1] - ai[i]) + 1; /* allocate some extra space */
190b05258aeSStefano Zampini         miss          = PETSC_TRUE;
191b05258aeSStefano Zampini       } else {
192421480d9SBarry Smith         rowlengths[i] = (ai[i + 1] - adiag[i]);
19359557b74SHong Zhang       }
194b05258aeSStefano Zampini     }
1955a2b941aSBarry Smith   } else PetscCall(MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(A, &rowlengths));
196bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
1979566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1989566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
1999566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSBAIJ));
2009566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, rowlengths));
201bd019fc1SStefano Zampini   } else B = *newmat;
20259557b74SHong Zhang 
203b05258aeSStefano Zampini   if (bs == 1 && !miss) {
204f4f49eeaSPierre Jolivet     b  = (Mat_SeqSBAIJ *)B->data;
205861ba921SHong Zhang     bi = b->i;
206861ba921SHong Zhang     bj = b->j;
207861ba921SHong Zhang     bv = b->a;
208861ba921SHong Zhang 
209861ba921SHong Zhang     bi[0] = 0;
21059557b74SHong Zhang     for (i = 0; i < m; i++) {
211421480d9SBarry Smith       aj = a->j + adiag[i];
212421480d9SBarry Smith       av = a->a + adiag[i];
213861ba921SHong Zhang       for (j = 0; j < rowlengths[i]; j++) {
2149371c9d4SSatish Balay         *bj = *aj;
2159371c9d4SSatish Balay         bj++;
2169371c9d4SSatish Balay         aj++;
2179371c9d4SSatish Balay         *bv = *av;
2189371c9d4SSatish Balay         bv++;
2199371c9d4SSatish Balay         av++;
220861ba921SHong Zhang       }
221861ba921SHong Zhang       bi[i + 1]  = bi[i] + rowlengths[i];
222861ba921SHong Zhang       b->ilen[i] = rowlengths[i];
22359557b74SHong Zhang     }
2249566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2259566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
226ae8d29abSPierre Jolivet   } else {
2279566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
228ae8d29abSPierre 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 */
229ae8d29abSPierre Jolivet     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
230ae8d29abSPierre Jolivet     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
2319566063dSJacob Faibussowitsch     PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B));
232ae8d29abSPierre Jolivet   }
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowlengths));
234ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
235ac530a7eSPierre Jolivet   else *newmat = B;
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23759557b74SHong Zhang }
23859557b74SHong Zhang 
239d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
240d71ae5a4SJacob Faibussowitsch {
241a0e1a404SHong Zhang   Mat           B;
242a0e1a404SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
243a0e1a404SHong Zhang   Mat_SeqBAIJ  *b;
244421480d9SBarry Smith   PetscInt     *ai = a->i, *aj = a->j, m = A->rmap->N, n = A->cmap->n, *bi, *bj, *browlengths, nz, *browstart, itmp;
245421480d9SBarry Smith   PetscInt      bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs;
246dd6ea824SBarry Smith   MatScalar    *av, *bv;
247eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
248b94d7dedSBarry Smith   const int aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
249eb1ec7c1SStefano Zampini #else
250eb1ec7c1SStefano Zampini   const int aconj = 0;
251eb1ec7c1SStefano Zampini #endif
252a0e1a404SHong Zhang 
253a0e1a404SHong Zhang   PetscFunctionBegin;
254a0e1a404SHong Zhang   /* compute browlengths of newmat */
2559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &browlengths, mbs, &browstart));
256421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) browlengths[i] = 0;
257421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
258a0e1a404SHong Zhang     nz = ai[i + 1] - ai[i];
259a0e1a404SHong Zhang     aj++;                               /* skip diagonal */
260421480d9SBarry Smith     for (PetscInt k = 1; k < nz; k++) { /* no. of lower triangular blocks */
2619371c9d4SSatish Balay       browlengths[*aj]++;
2629371c9d4SSatish Balay       aj++;
263a0e1a404SHong Zhang     }
264a0e1a404SHong Zhang     browlengths[i] += nz; /* no. of upper triangular blocks */
265a0e1a404SHong Zhang   }
266a0e1a404SHong Zhang 
267bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2689566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2699566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
2709566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQBAIJ));
2719566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetPreallocation(B, bs, 0, browlengths));
272bd019fc1SStefano Zampini   } else B = *newmat;
273a0e1a404SHong Zhang 
274f4f49eeaSPierre Jolivet   b  = (Mat_SeqBAIJ *)B->data;
275a0e1a404SHong Zhang   bi = b->i;
276a0e1a404SHong Zhang   bj = b->j;
277a0e1a404SHong Zhang   bv = b->a;
278a0e1a404SHong Zhang 
279a0e1a404SHong Zhang   /* set b->i */
280a0e1a404SHong Zhang   bi[0] = 0;
281421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
282a0e1a404SHong Zhang     b->ilen[i]   = browlengths[i];
283a0e1a404SHong Zhang     bi[i + 1]    = bi[i] + browlengths[i];
284a0e1a404SHong Zhang     browstart[i] = bi[i];
285a0e1a404SHong Zhang   }
286aed4548fSBarry 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);
287a0e1a404SHong Zhang 
288a0e1a404SHong Zhang   /* set b->j and b->a */
2899371c9d4SSatish Balay   aj = a->j;
2909371c9d4SSatish Balay   av = a->a;
291421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
292a0e1a404SHong Zhang     /* diagonal block */
2939371c9d4SSatish Balay     *(bj + browstart[i]) = *aj;
2949371c9d4SSatish Balay     aj++;
29526fbe8dcSKarl Rupp 
296a0e1a404SHong Zhang     itmp = bs2 * browstart[i];
297421480d9SBarry Smith     for (PetscInt k = 0; k < bs2; k++) {
2989371c9d4SSatish Balay       *(bv + itmp + k) = *av;
2999371c9d4SSatish Balay       av++;
300a0e1a404SHong Zhang     }
301a0e1a404SHong Zhang     browstart[i]++;
302a0e1a404SHong Zhang 
303a0e1a404SHong Zhang     nz = ai[i + 1] - ai[i] - 1;
304a0e1a404SHong Zhang     while (nz--) {
30574ee4d9fSHong Zhang       /* lower triangular blocks - transpose blocks of A */
30674ee4d9fSHong Zhang       *(bj + browstart[*aj]) = i; /* block col index */
30726fbe8dcSKarl Rupp 
30874ee4d9fSHong Zhang       itmp = bs2 * browstart[*aj]; /* row index */
309421480d9SBarry Smith       for (PetscInt col = 0; col < bs; col++) {
310421480d9SBarry Smith         PetscInt k = col;
311421480d9SBarry Smith         for (PetscInt row = 0; row < bs; row++) {
312eb1ec7c1SStefano Zampini           bv[itmp + col * bs + row] = aconj ? PetscConj(av[k]) : av[k];
313eb1ec7c1SStefano Zampini           k += bs;
31474ee4d9fSHong Zhang         }
315a0e1a404SHong Zhang       }
316a0e1a404SHong Zhang       browstart[*aj]++;
317a0e1a404SHong Zhang 
318a0e1a404SHong Zhang       /* upper triangular blocks */
3199371c9d4SSatish Balay       *(bj + browstart[i]) = *aj;
3209371c9d4SSatish Balay       aj++;
32126fbe8dcSKarl Rupp 
322a0e1a404SHong Zhang       itmp = bs2 * browstart[i];
323421480d9SBarry Smith       for (PetscInt k = 0; k < bs2; k++) bv[itmp + k] = av[k];
32474ee4d9fSHong Zhang       av += bs2;
325a0e1a404SHong Zhang       browstart[i]++;
326a0e1a404SHong Zhang     }
327a0e1a404SHong Zhang   }
3289566063dSJacob Faibussowitsch   PetscCall(PetscFree2(browlengths, browstart));
3299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
331a0e1a404SHong Zhang 
332ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
333ac530a7eSPierre Jolivet   else *newmat = B;
3343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
335a0e1a404SHong Zhang }
336be1d678aSKris Buschelman 
337d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
338d71ae5a4SJacob Faibussowitsch {
339a0e1a404SHong Zhang   Mat             B;
340a0e1a404SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data;
341a0e1a404SHong Zhang   Mat_SeqSBAIJ   *b;
342421480d9SBarry Smith   PetscInt       *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->n, k, *bi, *bj, *browlengths;
343421480d9SBarry Smith   PetscInt        bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs;
344dd6ea824SBarry Smith   MatScalar      *av, *bv;
345421480d9SBarry Smith   const PetscInt *adiag;
346a0e1a404SHong Zhang 
347a0e1a404SHong Zhang   PetscFunctionBegin;
348*b0c98d1dSPierre Jolivet   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
34908401ef6SPierre Jolivet   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square");
350421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &adiag, NULL));
3519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs, &browlengths));
352421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) browlengths[i] = ai[i + 1] - adiag[i];
353a0e1a404SHong Zhang 
354bd019fc1SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
3559566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
3569566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
3579566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSBAIJ));
3589566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, browlengths));
359bd019fc1SStefano Zampini   } else B = *newmat;
360a0e1a404SHong Zhang 
361f4f49eeaSPierre Jolivet   b  = (Mat_SeqSBAIJ *)B->data;
362a0e1a404SHong Zhang   bi = b->i;
363a0e1a404SHong Zhang   bj = b->j;
364a0e1a404SHong Zhang   bv = b->a;
365a0e1a404SHong Zhang 
366a0e1a404SHong Zhang   bi[0] = 0;
367421480d9SBarry Smith   for (PetscInt i = 0; i < mbs; i++) {
368421480d9SBarry Smith     aj = a->j + adiag[i];
369421480d9SBarry Smith     av = a->a + (adiag[i]) * bs2;
370421480d9SBarry Smith     for (PetscInt j = 0; j < browlengths[i]; j++) {
3719371c9d4SSatish Balay       *bj = *aj;
3729371c9d4SSatish Balay       bj++;
3739371c9d4SSatish Balay       aj++;
374a0e1a404SHong Zhang       for (k = 0; k < bs2; k++) {
3759371c9d4SSatish Balay         *bv = *av;
3769371c9d4SSatish Balay         bv++;
3779371c9d4SSatish Balay         av++;
378a0e1a404SHong Zhang       }
379a0e1a404SHong Zhang     }
380a0e1a404SHong Zhang     bi[i + 1]  = bi[i] + browlengths[i];
381a0e1a404SHong Zhang     b->ilen[i] = browlengths[i];
382a0e1a404SHong Zhang   }
3839566063dSJacob Faibussowitsch   PetscCall(PetscFree(browlengths));
3849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
386a0e1a404SHong Zhang 
387ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
388ac530a7eSPierre Jolivet   else *newmat = B;
3893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
390a0e1a404SHong Zhang }
391