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 52*f4f49eeaSPierre 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; 129f85a0629SBarry Smith PetscInt m, n, bs = PetscAbs(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) { 136f85a0629SBarry Smith const PetscInt *adiag = Aa->diag; 137f85a0629SBarry Smith 138f85a0629SBarry Smith PetscCall(PetscMalloc1(m, nnz)); 139f85a0629SBarry Smith for (PetscInt i = 0; i < m; i++) { 140f85a0629SBarry Smith if (adiag[i] == ai[i + 1]) { 141f85a0629SBarry Smith (*nnz)[i] = 0; 142f85a0629SBarry Smith for (PetscInt j = ai[i]; j < ai[i + 1]; j++) (*nnz)[i] += (aj[j] > i); 143f85a0629SBarry Smith } else (*nnz)[i] = ai[i + 1] - adiag[i]; 144f85a0629SBarry Smith } 145f85a0629SBarry Smith } else { 146f85a0629SBarry Smith PetscHSetIJ ht; 147f85a0629SBarry Smith PetscHashIJKey key; 148f85a0629SBarry Smith PetscBool missing; 149f85a0629SBarry Smith 150f85a0629SBarry Smith PetscCall(PetscHSetIJCreate(&ht)); 151f85a0629SBarry Smith PetscCall(PetscCalloc1(m / bs, nnz)); 152f85a0629SBarry Smith for (PetscInt i = 0; i < m; i++) { 153f85a0629SBarry Smith key.i = i / bs; 154f85a0629SBarry Smith for (PetscInt k = ai[i]; k < ai[i + 1]; k++) { 155f85a0629SBarry Smith key.j = aj[k] / bs; 156f85a0629SBarry Smith if (key.j < key.i) continue; 157f85a0629SBarry Smith PetscCall(PetscHSetIJQueryAdd(ht, key, &missing)); 158f85a0629SBarry Smith if (missing) (*nnz)[key.i]++; 1595a2b941aSBarry Smith } 1605a2b941aSBarry Smith } 161f85a0629SBarry Smith PetscCall(PetscHSetIJDestroy(&ht)); 1625a2b941aSBarry Smith } 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1645a2b941aSBarry Smith } 1655a2b941aSBarry Smith 166d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 167d71ae5a4SJacob Faibussowitsch { 168676c34cdSKris Buschelman Mat B; 16959557b74SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 170861ba921SHong Zhang Mat_SeqSBAIJ *b; 171ae8d29abSPierre Jolivet PetscInt *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->N, i, j, *bi, *bj, *rowlengths, bs = PetscAbs(A->rmap->bs); 172dd6ea824SBarry Smith MatScalar *av, *bv; 173b05258aeSStefano Zampini PetscBool miss = PETSC_FALSE; 17459557b74SHong Zhang 17559557b74SHong Zhang PetscFunctionBegin; 176b05258aeSStefano Zampini #if !defined(PETSC_USE_COMPLEX) 177b94d7dedSBarry Smith PetscCheck(A->symmetric == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); 178b05258aeSStefano Zampini #else 179b94d7dedSBarry 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)"); 180b05258aeSStefano Zampini #endif 18108401ef6SPierre Jolivet PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square"); 18259557b74SHong Zhang 1835a2b941aSBarry Smith if (bs == 1) { 1845a2b941aSBarry Smith PetscCall(PetscMalloc1(m, &rowlengths)); 1855a2b941aSBarry Smith for (i = 0; i < m; i++) { 1865a2b941aSBarry Smith if (a->diag[i] == ai[i + 1]) { /* missing diagonal */ 1875a2b941aSBarry Smith rowlengths[i] = (ai[i + 1] - ai[i]) + 1; /* allocate some extra space */ 188b05258aeSStefano Zampini miss = PETSC_TRUE; 189b05258aeSStefano Zampini } else { 1905a2b941aSBarry Smith rowlengths[i] = (ai[i + 1] - a->diag[i]); 19159557b74SHong Zhang } 192b05258aeSStefano Zampini } 1935a2b941aSBarry Smith } else PetscCall(MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(A, &rowlengths)); 194bd019fc1SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 1959566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 1969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, n, m, n)); 1979566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQSBAIJ)); 1989566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, rowlengths)); 199bd019fc1SStefano Zampini } else B = *newmat; 20059557b74SHong Zhang 201b05258aeSStefano Zampini if (bs == 1 && !miss) { 202*f4f49eeaSPierre Jolivet b = (Mat_SeqSBAIJ *)B->data; 203861ba921SHong Zhang bi = b->i; 204861ba921SHong Zhang bj = b->j; 205861ba921SHong Zhang bv = b->a; 206861ba921SHong Zhang 207861ba921SHong Zhang bi[0] = 0; 20859557b74SHong Zhang for (i = 0; i < m; i++) { 20959557b74SHong Zhang aj = a->j + a->diag[i]; 21059557b74SHong Zhang av = a->a + a->diag[i]; 211861ba921SHong Zhang for (j = 0; j < rowlengths[i]; j++) { 2129371c9d4SSatish Balay *bj = *aj; 2139371c9d4SSatish Balay bj++; 2149371c9d4SSatish Balay aj++; 2159371c9d4SSatish Balay *bv = *av; 2169371c9d4SSatish Balay bv++; 2179371c9d4SSatish Balay av++; 218861ba921SHong Zhang } 219861ba921SHong Zhang bi[i + 1] = bi[i] + rowlengths[i]; 220861ba921SHong Zhang b->ilen[i] = rowlengths[i]; 22159557b74SHong Zhang } 2229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 224ae8d29abSPierre Jolivet } else { 2259566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE)); 226ae8d29abSPierre 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 */ 227ae8d29abSPierre Jolivet /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before */ 228ae8d29abSPierre Jolivet /* MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below */ 2299566063dSJacob Faibussowitsch PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B)); 230ae8d29abSPierre Jolivet } 2319566063dSJacob Faibussowitsch PetscCall(PetscFree(rowlengths)); 232511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 2339566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 234ae8d29abSPierre Jolivet } else *newmat = B; 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23659557b74SHong Zhang } 23759557b74SHong Zhang 238d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 239d71ae5a4SJacob Faibussowitsch { 240a0e1a404SHong Zhang Mat B; 241a0e1a404SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data; 242a0e1a404SHong Zhang Mat_SeqBAIJ *b; 243d0f46423SBarry Smith PetscInt *ai = a->i, *aj = a->j, m = A->rmap->N, n = A->cmap->n, i, k, *bi, *bj, *browlengths, nz, *browstart, itmp; 24474ee4d9fSHong Zhang PetscInt bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs, col, row; 245dd6ea824SBarry Smith MatScalar *av, *bv; 246eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 247b94d7dedSBarry Smith const int aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0; 248eb1ec7c1SStefano Zampini #else 249eb1ec7c1SStefano Zampini const int aconj = 0; 250eb1ec7c1SStefano Zampini #endif 251a0e1a404SHong Zhang 252a0e1a404SHong Zhang PetscFunctionBegin; 253a0e1a404SHong Zhang /* compute browlengths of newmat */ 2549566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mbs, &browlengths, mbs, &browstart)); 255a0e1a404SHong Zhang for (i = 0; i < mbs; i++) browlengths[i] = 0; 256a0e1a404SHong Zhang for (i = 0; i < mbs; i++) { 257a0e1a404SHong Zhang nz = ai[i + 1] - ai[i]; 258a0e1a404SHong Zhang aj++; /* skip diagonal */ 259a0e1a404SHong Zhang for (k = 1; k < nz; k++) { /* no. of lower triangular blocks */ 2609371c9d4SSatish Balay browlengths[*aj]++; 2619371c9d4SSatish Balay aj++; 262a0e1a404SHong Zhang } 263a0e1a404SHong Zhang browlengths[i] += nz; /* no. of upper triangular blocks */ 264a0e1a404SHong Zhang } 265a0e1a404SHong Zhang 266bd019fc1SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 2679566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 2689566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, n, m, n)); 2699566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQBAIJ)); 2709566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(B, bs, 0, browlengths)); 271bd019fc1SStefano Zampini } else B = *newmat; 272a0e1a404SHong Zhang 273*f4f49eeaSPierre Jolivet b = (Mat_SeqBAIJ *)B->data; 274a0e1a404SHong Zhang bi = b->i; 275a0e1a404SHong Zhang bj = b->j; 276a0e1a404SHong Zhang bv = b->a; 277a0e1a404SHong Zhang 278a0e1a404SHong Zhang /* set b->i */ 279a0e1a404SHong Zhang bi[0] = 0; 280a0e1a404SHong Zhang for (i = 0; i < mbs; i++) { 281a0e1a404SHong Zhang b->ilen[i] = browlengths[i]; 282a0e1a404SHong Zhang bi[i + 1] = bi[i] + browlengths[i]; 283a0e1a404SHong Zhang browstart[i] = bi[i]; 284a0e1a404SHong Zhang } 285aed4548fSBarry 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); 286a0e1a404SHong Zhang 287a0e1a404SHong Zhang /* set b->j and b->a */ 2889371c9d4SSatish Balay aj = a->j; 2899371c9d4SSatish Balay av = a->a; 290a0e1a404SHong Zhang for (i = 0; i < mbs; i++) { 291a0e1a404SHong Zhang /* diagonal block */ 2929371c9d4SSatish Balay *(bj + browstart[i]) = *aj; 2939371c9d4SSatish Balay aj++; 29426fbe8dcSKarl Rupp 295a0e1a404SHong Zhang itmp = bs2 * browstart[i]; 296a0e1a404SHong Zhang for (k = 0; k < bs2; k++) { 2979371c9d4SSatish Balay *(bv + itmp + k) = *av; 2989371c9d4SSatish Balay av++; 299a0e1a404SHong Zhang } 300a0e1a404SHong Zhang browstart[i]++; 301a0e1a404SHong Zhang 302a0e1a404SHong Zhang nz = ai[i + 1] - ai[i] - 1; 303a0e1a404SHong Zhang while (nz--) { 30474ee4d9fSHong Zhang /* lower triangular blocks - transpose blocks of A */ 30574ee4d9fSHong Zhang *(bj + browstart[*aj]) = i; /* block col index */ 30626fbe8dcSKarl Rupp 30774ee4d9fSHong Zhang itmp = bs2 * browstart[*aj]; /* row index */ 30874ee4d9fSHong Zhang for (col = 0; col < bs; col++) { 30974ee4d9fSHong Zhang k = col; 31074ee4d9fSHong Zhang for (row = 0; row < bs; row++) { 311eb1ec7c1SStefano Zampini bv[itmp + col * bs + row] = aconj ? PetscConj(av[k]) : av[k]; 312eb1ec7c1SStefano Zampini k += bs; 31374ee4d9fSHong Zhang } 314a0e1a404SHong Zhang } 315a0e1a404SHong Zhang browstart[*aj]++; 316a0e1a404SHong Zhang 317a0e1a404SHong Zhang /* upper triangular blocks */ 3189371c9d4SSatish Balay *(bj + browstart[i]) = *aj; 3199371c9d4SSatish Balay aj++; 32026fbe8dcSKarl Rupp 321a0e1a404SHong Zhang itmp = bs2 * browstart[i]; 322ad540459SPierre Jolivet for (k = 0; k < bs2; k++) bv[itmp + k] = av[k]; 32374ee4d9fSHong Zhang av += bs2; 324a0e1a404SHong Zhang browstart[i]++; 325a0e1a404SHong Zhang } 326a0e1a404SHong Zhang } 3279566063dSJacob Faibussowitsch PetscCall(PetscFree2(browlengths, browstart)); 3289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 3299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 330a0e1a404SHong Zhang 331511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 3329566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 333ae8d29abSPierre 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; 342d0f46423SBarry Smith PetscInt *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->n, i, j, k, *bi, *bj, *browlengths; 343d0f46423SBarry Smith PetscInt bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs, dd; 344dd6ea824SBarry Smith MatScalar *av, *bv; 345ace3abfcSBarry Smith PetscBool flg; 346a0e1a404SHong Zhang 347a0e1a404SHong Zhang PetscFunctionBegin; 34828b400f6SJacob Faibussowitsch PetscCheck(A->symmetric, 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"); 3509566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal_SeqBAIJ(A, &flg, &dd)); /* check for missing diagonals, then mark diag */ 35128b400f6SJacob Faibussowitsch PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal %" PetscInt_FMT, dd); 352a0e1a404SHong Zhang 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs, &browlengths)); 354ad540459SPierre Jolivet for (i = 0; i < mbs; i++) browlengths[i] = ai[i + 1] - a->diag[i]; 355a0e1a404SHong Zhang 356bd019fc1SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 3579566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 3589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, n, m, n)); 3599566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQSBAIJ)); 3609566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, browlengths)); 361bd019fc1SStefano Zampini } else B = *newmat; 362a0e1a404SHong Zhang 363*f4f49eeaSPierre Jolivet b = (Mat_SeqSBAIJ *)B->data; 364a0e1a404SHong Zhang bi = b->i; 365a0e1a404SHong Zhang bj = b->j; 366a0e1a404SHong Zhang bv = b->a; 367a0e1a404SHong Zhang 368a0e1a404SHong Zhang bi[0] = 0; 369a0e1a404SHong Zhang for (i = 0; i < mbs; i++) { 370a0e1a404SHong Zhang aj = a->j + a->diag[i]; 371a0e1a404SHong Zhang av = a->a + (a->diag[i]) * bs2; 372a0e1a404SHong Zhang for (j = 0; j < browlengths[i]; j++) { 3739371c9d4SSatish Balay *bj = *aj; 3749371c9d4SSatish Balay bj++; 3759371c9d4SSatish Balay aj++; 376a0e1a404SHong Zhang for (k = 0; k < bs2; k++) { 3779371c9d4SSatish Balay *bv = *av; 3789371c9d4SSatish Balay bv++; 3799371c9d4SSatish Balay av++; 380a0e1a404SHong Zhang } 381a0e1a404SHong Zhang } 382a0e1a404SHong Zhang bi[i + 1] = bi[i] + browlengths[i]; 383a0e1a404SHong Zhang b->ilen[i] = browlengths[i]; 384a0e1a404SHong Zhang } 3859566063dSJacob Faibussowitsch PetscCall(PetscFree(browlengths)); 3869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 3879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 388a0e1a404SHong Zhang 389511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 3909566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 391ae8d29abSPierre Jolivet } else *newmat = B; 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 393a0e1a404SHong Zhang } 394