1be1d678aSKris Buschelman #define PETSCMAT_DLL 259557b74SHong Zhang 359557b74SHong Zhang #include "src/mat/impls/aij/seq/aij.h" 4595208bcSHong Zhang #include "src/mat/impls/baij/seq/baij.h" 5861ba921SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 659557b74SHong Zhang 759557b74SHong Zhang EXTERN_C_BEGIN 859557b74SHong Zhang #undef __FUNCT__ 94e5e7fe4SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqAIJ" 10f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 114e5e7fe4SHong Zhang { 124e5e7fe4SHong Zhang Mat B; 134e5e7fe4SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 144e5e7fe4SHong Zhang Mat_SeqAIJ *b; 15dfbe8321SBarry Smith PetscErrorCode ierr; 16899cda47SBarry 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; 17899cda47SBarry Smith PetscInt bs=A->rmap.bs,bs2=bs*bs,mbs=A->rmap.N/bs; 184e5e7fe4SHong Zhang PetscScalar *av,*bv; 194e5e7fe4SHong Zhang 204e5e7fe4SHong Zhang PetscFunctionBegin; 214e5e7fe4SHong Zhang /* compute rowlengths of newmat */ 2213f74950SBarry Smith ierr = PetscMalloc((2*m+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 23a7a3a9ebSHong Zhang rowstart = rowlengths + m; 24a7a3a9ebSHong Zhang 25a7a3a9ebSHong Zhang for (i=0; i<mbs; i++) rowlengths[i*bs] = 0; 264e5e7fe4SHong Zhang aj = a->j; 27a7a3a9ebSHong Zhang k = 0; 28a7a3a9ebSHong Zhang for (i=0; i<mbs; i++) { 294e5e7fe4SHong Zhang nz = ai[i+1] - ai[i]; 30a7a3a9ebSHong Zhang aj++; /* skip diagonal */ 31a7a3a9ebSHong Zhang for (j=1; j<nz; j++) { /* no. of lower triangular blocks */ 32a7a3a9ebSHong Zhang rowlengths[(*aj)*bs]++; aj++; 33a7a3a9ebSHong Zhang } 34a7a3a9ebSHong Zhang rowlengths[k] += nz; /* no. of upper triangular blocks */ 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 43f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 44f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 45f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 46f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 47*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 48*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROWS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 49*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_COLUMNS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 50899cda47SBarry Smith B->rmap.bs = A->rmap.bs; 514e5e7fe4SHong Zhang 524e5e7fe4SHong Zhang 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 */ 58a7a3a9ebSHong Zhang bi[0] = 0; rowstart[0] = 0; 59a7a3a9ebSHong Zhang for (i=0; i<mbs; i++){ 60a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ 61a7a3a9ebSHong Zhang b->ilen[i*bs+j] = rowlengths[i*bs]; 62a7a3a9ebSHong Zhang rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs]; 634e5e7fe4SHong Zhang } 64a7a3a9ebSHong Zhang bi[i+1] = bi[i] + rowlengths[i*bs]/bs; 65a7a3a9ebSHong Zhang } 6677431f27SBarry Smith if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-mbs: %D\n",bi[mbs],2*a->nz - mbs); 674e5e7fe4SHong Zhang 684e5e7fe4SHong Zhang /* set b->j and b->a */ 694e5e7fe4SHong Zhang aj = a->j; av = a->a; 70a7a3a9ebSHong Zhang for (i=0; i<mbs; i++) { 71a7a3a9ebSHong Zhang /* diagonal block */ 72a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row i*bs+j */ 73a7a3a9ebSHong Zhang itmp = i*bs+j; 74a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col i*bs+k */ 75a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = (*aj)*bs+k; 76a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 77a7a3a9ebSHong Zhang rowstart[itmp]++; 78a7a3a9ebSHong Zhang } 79a7a3a9ebSHong Zhang } 80a7a3a9ebSHong Zhang aj++; av += bs2; 81a7a3a9ebSHong Zhang 824e5e7fe4SHong Zhang nz = ai[i+1] - ai[i] -1; 834e5e7fe4SHong Zhang while (nz--){ 84a7a3a9ebSHong Zhang /* lower triangular blocks */ 85a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row (*aj)*bs+j */ 86a7a3a9ebSHong Zhang itmp = (*aj)*bs+j; 87a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col i*bs+k */ 88a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = i*bs+k; 89a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 90a7a3a9ebSHong Zhang rowstart[itmp]++; 91a7a3a9ebSHong Zhang } 92a7a3a9ebSHong Zhang } 93a7a3a9ebSHong Zhang /* upper triangular blocks */ 94a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row i*bs+j */ 95a7a3a9ebSHong Zhang itmp = i*bs+j; 96a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col (*aj)*bs+k */ 97a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = (*aj)*bs+k; 98a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 99a7a3a9ebSHong Zhang rowstart[itmp]++; 100a7a3a9ebSHong Zhang } 101a7a3a9ebSHong Zhang } 102a7a3a9ebSHong Zhang aj++; av += bs2; 1034e5e7fe4SHong Zhang } 1044e5e7fe4SHong Zhang } 1054e5e7fe4SHong Zhang ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1064e5e7fe4SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1074e5e7fe4SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1084e5e7fe4SHong Zhang 109ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 1108ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 111c3d102feSKris Buschelman } else { 1124e5e7fe4SHong Zhang *newmat = B; 113c3d102feSKris Buschelman } 1144e5e7fe4SHong Zhang PetscFunctionReturn(0); 1154e5e7fe4SHong Zhang } 116be1d678aSKris Buschelman EXTERN_C_END 117be1d678aSKris Buschelman 118be1d678aSKris Buschelman EXTERN_C_BEGIN 1194e5e7fe4SHong Zhang #undef __FUNCT__ 12059557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ" 121f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) { 122676c34cdSKris Buschelman Mat B; 12359557b74SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 124861ba921SHong Zhang Mat_SeqSBAIJ *b; 125dfbe8321SBarry Smith PetscErrorCode ierr; 126899cda47SBarry Smith PetscInt *ai=a->i,*aj,m=A->rmap.N,n=A->cmap.N,i,j,*bi,*bj,*rowlengths; 127861ba921SHong Zhang PetscScalar *av,*bv; 12859557b74SHong Zhang 12959557b74SHong Zhang PetscFunctionBegin; 1302d9a3abdSHong Zhang if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 13159557b74SHong Zhang 13213f74950SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 13359557b74SHong Zhang for (i=0; i<m; i++) { 13459557b74SHong Zhang rowlengths[i] = ai[i+1] - a->diag[i]; 13559557b74SHong Zhang } 136f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 137f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 138f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 139ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr); 14059557b74SHong Zhang 141*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 142*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROWS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 143*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_COLUMNS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 14459557b74SHong Zhang 145676c34cdSKris Buschelman b = (Mat_SeqSBAIJ*)(B->data); 146861ba921SHong Zhang bi = b->i; 147861ba921SHong Zhang bj = b->j; 148861ba921SHong Zhang bv = b->a; 149861ba921SHong Zhang 150861ba921SHong Zhang bi[0] = 0; 15159557b74SHong Zhang for (i=0; i<m; i++) { 15259557b74SHong Zhang aj = a->j + a->diag[i]; 15359557b74SHong Zhang av = a->a + a->diag[i]; 154861ba921SHong Zhang for (j=0; j<rowlengths[i]; j++){ 155861ba921SHong Zhang *bj = *aj; bj++; aj++; 156861ba921SHong Zhang *bv = *av; bv++; av++; 157861ba921SHong Zhang } 158861ba921SHong Zhang bi[i+1] = bi[i] + rowlengths[i]; 159861ba921SHong Zhang b->ilen[i] = rowlengths[i]; 16059557b74SHong Zhang } 16159557b74SHong Zhang 16259557b74SHong Zhang ierr = PetscFree(rowlengths);CHKERRQ(ierr); 163676c34cdSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164676c34cdSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165676c34cdSKris Buschelman 166ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 1678ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 168c3d102feSKris Buschelman } else { 169676c34cdSKris Buschelman *newmat = B; 170c3d102feSKris Buschelman } 17159557b74SHong Zhang PetscFunctionReturn(0); 17259557b74SHong Zhang } 17359557b74SHong Zhang EXTERN_C_END 17459557b74SHong Zhang 175a0e1a404SHong Zhang EXTERN_C_BEGIN 176a0e1a404SHong Zhang #undef __FUNCT__ 177a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ" 178f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 179a0e1a404SHong Zhang { 180a0e1a404SHong Zhang Mat B; 181a0e1a404SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 182a0e1a404SHong Zhang Mat_SeqBAIJ *b; 183dfbe8321SBarry Smith PetscErrorCode ierr; 184899cda47SBarry Smith PetscInt *ai=a->i,*aj=a->j,m=A->rmap.N,n=A->cmap.n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp; 185899cda47SBarry Smith PetscInt bs=A->rmap.bs,bs2=bs*bs,mbs=m/bs; 186a0e1a404SHong Zhang PetscScalar *av,*bv; 187a0e1a404SHong Zhang 188a0e1a404SHong Zhang PetscFunctionBegin; 189a0e1a404SHong Zhang /* compute browlengths of newmat */ 19013f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); 191a0e1a404SHong Zhang browstart = browlengths + mbs; 192a0e1a404SHong Zhang for (i=0; i<mbs; i++) browlengths[i] = 0; 193a0e1a404SHong Zhang aj = a->j; 194a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 195a0e1a404SHong Zhang nz = ai[i+1] - ai[i]; 196a0e1a404SHong Zhang aj++; /* skip diagonal */ 197a0e1a404SHong Zhang for (k=1; k<nz; k++) { /* no. of lower triangular blocks */ 198a0e1a404SHong Zhang browlengths[*aj]++; aj++; 199a0e1a404SHong Zhang } 200a0e1a404SHong Zhang browlengths[i] += nz; /* no. of upper triangular blocks */ 201a0e1a404SHong Zhang } 202a0e1a404SHong Zhang 203f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 204f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 205f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 206f204ca49SKris Buschelman ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 207*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 208*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROWS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 209*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_COLUMNS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 210a0e1a404SHong Zhang 211a0e1a404SHong Zhang b = (Mat_SeqBAIJ*)(B->data); 212a0e1a404SHong Zhang bi = b->i; 213a0e1a404SHong Zhang bj = b->j; 214a0e1a404SHong Zhang bv = b->a; 215a0e1a404SHong Zhang 216a0e1a404SHong Zhang /* set b->i */ 217a0e1a404SHong Zhang bi[0] = 0; 218a0e1a404SHong Zhang for (i=0; i<mbs; i++){ 219a0e1a404SHong Zhang b->ilen[i] = browlengths[i]; 220a0e1a404SHong Zhang bi[i+1] = bi[i] + browlengths[i]; 221a0e1a404SHong Zhang browstart[i] = bi[i]; 222a0e1a404SHong Zhang } 22377431f27SBarry Smith if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs); 224a0e1a404SHong Zhang 225a0e1a404SHong Zhang /* set b->j and b->a */ 226a0e1a404SHong Zhang aj = a->j; av = a->a; 227a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 228a0e1a404SHong Zhang /* diagonal block */ 229a0e1a404SHong Zhang *(bj + browstart[i]) = *aj; aj++; 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--){ 238a0e1a404SHong Zhang /* lower triangular blocks */ 239a0e1a404SHong Zhang *(bj + browstart[*aj]) = i; 240a0e1a404SHong Zhang itmp = bs2*browstart[*aj]; 241a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 242a0e1a404SHong Zhang *(bv + itmp + k) = *(av + k); 243a0e1a404SHong Zhang } 244a0e1a404SHong Zhang browstart[*aj]++; 245a0e1a404SHong Zhang 246a0e1a404SHong Zhang /* upper triangular blocks */ 247a0e1a404SHong Zhang *(bj + browstart[i]) = *aj; aj++; 248a0e1a404SHong Zhang itmp = bs2*browstart[i]; 249a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 250a0e1a404SHong Zhang *(bv + itmp + k) = *av; av++; 251a0e1a404SHong Zhang } 252a0e1a404SHong Zhang browstart[i]++; 253a0e1a404SHong Zhang } 254a0e1a404SHong Zhang } 255a0e1a404SHong Zhang ierr = PetscFree(browlengths);CHKERRQ(ierr); 256a0e1a404SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 257a0e1a404SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 258a0e1a404SHong Zhang 259ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 2608ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 261c3d102feSKris Buschelman } else { 262a0e1a404SHong Zhang *newmat = B; 263c3d102feSKris Buschelman } 264a0e1a404SHong Zhang PetscFunctionReturn(0); 265a0e1a404SHong Zhang } 266be1d678aSKris Buschelman EXTERN_C_END 267be1d678aSKris Buschelman 268be1d678aSKris Buschelman EXTERN_C_BEGIN 269a0e1a404SHong Zhang #undef __FUNCT__ 270a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ" 271f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 272a0e1a404SHong Zhang { 273a0e1a404SHong Zhang Mat B; 274a0e1a404SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 275a0e1a404SHong Zhang Mat_SeqSBAIJ *b; 276dfbe8321SBarry Smith PetscErrorCode ierr; 277899cda47SBarry Smith PetscInt *ai=a->i,*aj,m=A->rmap.N,n=A->cmap.n,i,j,k,*bi,*bj,*browlengths; 278899cda47SBarry Smith PetscInt bs=A->rmap.bs,bs2=bs*bs,mbs=m/bs; 279a0e1a404SHong Zhang PetscScalar *av,*bv; 280a0e1a404SHong Zhang 281a0e1a404SHong Zhang PetscFunctionBegin; 282a0e1a404SHong Zhang if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 283595208bcSHong Zhang ierr = MatMissingDiagonal_SeqBAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ 284a0e1a404SHong Zhang 28513f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); 286a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 287a0e1a404SHong Zhang browlengths[i] = ai[i+1] - a->diag[i]; 288a0e1a404SHong Zhang } 289a0e1a404SHong Zhang 290f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 291f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 292f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 293ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr); 294*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 295*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROWS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 296*4e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_COLUMNS_SORTED,PETSC_TRUE);CHKERRQ(ierr); 297a0e1a404SHong Zhang 298a0e1a404SHong Zhang b = (Mat_SeqSBAIJ*)(B->data); 299a0e1a404SHong Zhang bi = b->i; 300a0e1a404SHong Zhang bj = b->j; 301a0e1a404SHong Zhang bv = b->a; 302a0e1a404SHong Zhang 303a0e1a404SHong Zhang bi[0] = 0; 304a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 305a0e1a404SHong Zhang aj = a->j + a->diag[i]; 306a0e1a404SHong Zhang av = a->a + (a->diag[i])*bs2; 307a0e1a404SHong Zhang for (j=0; j<browlengths[i]; j++){ 308a0e1a404SHong Zhang *bj = *aj; bj++; aj++; 309a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 310a0e1a404SHong Zhang *bv = *av; bv++; av++; 311a0e1a404SHong Zhang } 312a0e1a404SHong Zhang } 313a0e1a404SHong Zhang bi[i+1] = bi[i] + browlengths[i]; 314a0e1a404SHong Zhang b->ilen[i] = browlengths[i]; 315a0e1a404SHong Zhang } 316a0e1a404SHong Zhang ierr = PetscFree(browlengths);CHKERRQ(ierr); 317a0e1a404SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 318a0e1a404SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 319a0e1a404SHong Zhang 320ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 3218ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 322c3d102feSKris Buschelman } else { 323a0e1a404SHong Zhang *newmat = B; 324c3d102feSKris Buschelman } 325a0e1a404SHong Zhang 326a0e1a404SHong Zhang PetscFunctionReturn(0); 327a0e1a404SHong Zhang } 328a0e1a404SHong Zhang EXTERN_C_END 329