1be1d678aSKris Buschelman #define PETSCMAT_DLL 259557b74SHong Zhang 37c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 47c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" 57c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 659557b74SHong Zhang 759557b74SHong Zhang EXTERN_C_BEGIN 859557b74SHong Zhang #undef __FUNCT__ 94124e67cSShri Abhyankar #define __FUNCT__ "MatConvert_SeqSBAIJ_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; 16d0f46423SBarry 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; 1701be0148SBarry Smith PetscInt bs=A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0; 18dd6ea824SBarry Smith MatScalar *av,*bv; 194e5e7fe4SHong Zhang 204e5e7fe4SHong Zhang PetscFunctionBegin; 214e5e7fe4SHong Zhang /* compute rowlengths of newmat */ 2274ed9c26SBarry Smith ierr = PetscMalloc2(m,PetscInt,&rowlengths,m+1,PetscInt,&rowstart);CHKERRQ(ierr); 23a7a3a9ebSHong Zhang 24a7a3a9ebSHong Zhang for (i=0; i<mbs; i++) rowlengths[i*bs] = 0; 254e5e7fe4SHong Zhang aj = a->j; 26a7a3a9ebSHong Zhang k = 0; 27a7a3a9ebSHong Zhang for (i=0; i<mbs; i++) { 284e5e7fe4SHong Zhang nz = ai[i+1] - ai[i]; 2901be0148SBarry Smith if (nz) { 3001be0148SBarry Smith rowlengths[k] += nz; /* no. of upper triangular blocks */ 3101be0148SBarry Smith if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */ 3201be0148SBarry Smith for (j=0; j<nz; j++) { /* no. of lower triangular blocks */ 33a7a3a9ebSHong Zhang rowlengths[(*aj)*bs]++; aj++; 34a7a3a9ebSHong Zhang } 3501be0148SBarry Smith } 36a7a3a9ebSHong Zhang rowlengths[k] *= bs; 37a7a3a9ebSHong Zhang for (j=1; j<bs; j++) { 38a7a3a9ebSHong Zhang rowlengths[k+j] = rowlengths[k]; 39a7a3a9ebSHong Zhang } 40a7a3a9ebSHong Zhang k += bs; 41a7a3a9ebSHong Zhang /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */ 424e5e7fe4SHong Zhang } 434e5e7fe4SHong Zhang 447adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 45f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 46f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 47f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 484e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 49d0f46423SBarry Smith B->rmap->bs = A->rmap->bs; 504e5e7fe4SHong Zhang 514e5e7fe4SHong Zhang b = (Mat_SeqAIJ*)(B->data); 524e5e7fe4SHong Zhang bi = b->i; 534e5e7fe4SHong Zhang bj = b->j; 544e5e7fe4SHong Zhang bv = b->a; 554e5e7fe4SHong Zhang 564e5e7fe4SHong Zhang /* set b->i */ 57a7a3a9ebSHong Zhang bi[0] = 0; rowstart[0] = 0; 58a7a3a9ebSHong Zhang for (i=0; i<mbs; i++){ 59a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ 60a7a3a9ebSHong Zhang b->ilen[i*bs+j] = rowlengths[i*bs]; 61a7a3a9ebSHong Zhang rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs]; 624e5e7fe4SHong Zhang } 63a7a3a9ebSHong Zhang bi[i+1] = bi[i] + rowlengths[i*bs]/bs; 64a7a3a9ebSHong Zhang } 6501be0148SBarry Smith if (bi[mbs] != 2*a->nz - diagcnt) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-diagcnt: %D\n",bi[mbs],2*a->nz - diagcnt); 664e5e7fe4SHong Zhang 674e5e7fe4SHong Zhang /* set b->j and b->a */ 684e5e7fe4SHong Zhang aj = a->j; av = a->a; 69a7a3a9ebSHong Zhang for (i=0; i<mbs; i++) { 7001be0148SBarry Smith nz = ai[i+1] - ai[i]; 71a7a3a9ebSHong Zhang /* diagonal block */ 7201be0148SBarry Smith if (nz && *aj == i) { 7301be0148SBarry Smith nz--; 74a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row i*bs+j */ 75a7a3a9ebSHong Zhang itmp = i*bs+j; 76a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col i*bs+k */ 77a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = (*aj)*bs+k; 78a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 79a7a3a9ebSHong Zhang rowstart[itmp]++; 80a7a3a9ebSHong Zhang } 81a7a3a9ebSHong Zhang } 82a7a3a9ebSHong Zhang aj++; 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; 91*3bededecSBarry Smith *(bv + rowstart[itmp]) = *(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 } 104a7a3a9ebSHong Zhang aj++; av += bs2; 1054e5e7fe4SHong Zhang } 1064e5e7fe4SHong Zhang } 10774ed9c26SBarry Smith ierr = PetscFree2(rowlengths,rowstart);CHKERRQ(ierr); 1084e5e7fe4SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1094e5e7fe4SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1104e5e7fe4SHong Zhang 111ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 1128ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 113c3d102feSKris Buschelman } else { 1144e5e7fe4SHong Zhang *newmat = B; 115c3d102feSKris Buschelman } 1164e5e7fe4SHong Zhang PetscFunctionReturn(0); 1174e5e7fe4SHong Zhang } 118be1d678aSKris Buschelman EXTERN_C_END 119be1d678aSKris Buschelman 120be1d678aSKris Buschelman EXTERN_C_BEGIN 1214e5e7fe4SHong Zhang #undef __FUNCT__ 12259557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ" 1238cf70c4bSSatish Balay PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) { 124676c34cdSKris Buschelman Mat B; 12559557b74SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 126861ba921SHong Zhang Mat_SeqSBAIJ *b; 127dfbe8321SBarry Smith PetscErrorCode ierr; 128d0f46423SBarry Smith PetscInt *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths; 129dd6ea824SBarry Smith MatScalar *av,*bv; 13059557b74SHong Zhang 13159557b74SHong Zhang PetscFunctionBegin; 132e32f2f54SBarry Smith if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); 13359557b74SHong Zhang 13413f74950SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 13559557b74SHong Zhang for (i=0; i<m; i++) { 13659557b74SHong Zhang rowlengths[i] = ai[i+1] - a->diag[i]; 13759557b74SHong Zhang } 1387adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 139f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 140f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 141ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr); 14259557b74SHong Zhang 1434e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,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); 165eeffb40dSHong Zhang if (A->hermitian){ 166eeffb40dSHong Zhang ierr = MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); 167eeffb40dSHong Zhang } 168676c34cdSKris Buschelman 169ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 1708ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 171c3d102feSKris Buschelman } else { 172676c34cdSKris Buschelman *newmat = B; 173c3d102feSKris Buschelman } 17459557b74SHong Zhang PetscFunctionReturn(0); 17559557b74SHong Zhang } 17659557b74SHong Zhang EXTERN_C_END 17759557b74SHong Zhang 178a0e1a404SHong Zhang EXTERN_C_BEGIN 179a0e1a404SHong Zhang #undef __FUNCT__ 1807af9e4fdSHong Zhang #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqBAIJ" 181f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 182a0e1a404SHong Zhang { 183a0e1a404SHong Zhang Mat B; 184a0e1a404SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 185a0e1a404SHong Zhang Mat_SeqBAIJ *b; 186dfbe8321SBarry Smith PetscErrorCode ierr; 187d0f46423SBarry Smith PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp; 188d0f46423SBarry Smith PetscInt bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs; 189dd6ea824SBarry Smith MatScalar *av,*bv; 190a0e1a404SHong Zhang 191a0e1a404SHong Zhang PetscFunctionBegin; 192a0e1a404SHong Zhang /* compute browlengths of newmat */ 19374ed9c26SBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);CHKERRQ(ierr); 194a0e1a404SHong Zhang for (i=0; i<mbs; i++) browlengths[i] = 0; 195a0e1a404SHong Zhang aj = a->j; 196a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 197a0e1a404SHong Zhang nz = ai[i+1] - ai[i]; 198a0e1a404SHong Zhang aj++; /* skip diagonal */ 199a0e1a404SHong Zhang for (k=1; k<nz; k++) { /* no. of lower triangular blocks */ 200a0e1a404SHong Zhang browlengths[*aj]++; aj++; 201a0e1a404SHong Zhang } 202a0e1a404SHong Zhang browlengths[i] += nz; /* no. of upper triangular blocks */ 203a0e1a404SHong Zhang } 204a0e1a404SHong Zhang 2057adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 206f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 207f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 208f204ca49SKris Buschelman ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 2094e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,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 } 223e32f2f54SBarry Smith if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_COMM_SELF,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 } 25574ed9c26SBarry Smith ierr = PetscFree2(browlengths,browstart);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; 277d0f46423SBarry Smith PetscInt *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths; 278d0f46423SBarry Smith PetscInt bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd; 279dd6ea824SBarry Smith MatScalar *av,*bv; 280ace3abfcSBarry Smith PetscBool flg; 281a0e1a404SHong Zhang 282a0e1a404SHong Zhang PetscFunctionBegin; 283e32f2f54SBarry Smith if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); 2842af78befSBarry Smith ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ 285e32f2f54SBarry Smith if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd); 286a0e1a404SHong Zhang 28713f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); 288a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 289a0e1a404SHong Zhang browlengths[i] = ai[i+1] - a->diag[i]; 290a0e1a404SHong Zhang } 291a0e1a404SHong Zhang 2927adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 293f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 294f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 295ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr); 2964e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,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