159557b74SHong Zhang 259557b74SHong Zhang #include "src/mat/impls/aij/seq/aij.h" 3595208bcSHong Zhang #include "src/mat/impls/baij/seq/baij.h" 4861ba921SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 559557b74SHong Zhang 659557b74SHong Zhang EXTERN_C_BEGIN 759557b74SHong Zhang #undef __FUNCT__ 84e5e7fe4SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqAIJ" 9ceb03754SKris Buschelman PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) 104e5e7fe4SHong Zhang { 114e5e7fe4SHong Zhang Mat B; 124e5e7fe4SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 134e5e7fe4SHong Zhang Mat_SeqAIJ *b; 14dfbe8321SBarry Smith PetscErrorCode ierr; 1513f74950SBarry Smith PetscInt *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,j,k,*bi,*bj,*rowlengths,nz,*rowstart,itmp; 1613f74950SBarry Smith PetscInt bs=A->bs,bs2=bs*bs,mbs=A->m/bs; 174e5e7fe4SHong Zhang PetscScalar *av,*bv; 184e5e7fe4SHong Zhang 194e5e7fe4SHong Zhang PetscFunctionBegin; 20a7a3a9ebSHong Zhang 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 43f204ca49SKris Buschelman ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 44f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 45f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 46a7a3a9ebSHong Zhang ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr); 474e5e7fe4SHong Zhang ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 484e5e7fe4SHong Zhang ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 49521d7252SBarry Smith B->bs = A->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 } 6577431f27SBarry 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); 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++) { 70a7a3a9ebSHong Zhang /* diagonal block */ 71a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row i*bs+j */ 72a7a3a9ebSHong Zhang itmp = i*bs+j; 73a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col i*bs+k */ 74a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = (*aj)*bs+k; 75a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 76a7a3a9ebSHong Zhang rowstart[itmp]++; 77a7a3a9ebSHong Zhang } 78a7a3a9ebSHong Zhang } 79a7a3a9ebSHong Zhang aj++; av += bs2; 80a7a3a9ebSHong Zhang 814e5e7fe4SHong Zhang nz = ai[i+1] - ai[i] -1; 824e5e7fe4SHong Zhang while (nz--){ 83a7a3a9ebSHong Zhang /* lower triangular blocks */ 84a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row (*aj)*bs+j */ 85a7a3a9ebSHong Zhang itmp = (*aj)*bs+j; 86a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col i*bs+k */ 87a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = i*bs+k; 88a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 89a7a3a9ebSHong Zhang rowstart[itmp]++; 90a7a3a9ebSHong Zhang } 91a7a3a9ebSHong Zhang } 92a7a3a9ebSHong Zhang /* upper triangular blocks */ 93a7a3a9ebSHong Zhang for (j=0; j<bs; j++){ /* row i*bs+j */ 94a7a3a9ebSHong Zhang itmp = i*bs+j; 95a7a3a9ebSHong Zhang for (k=0; k<bs; k++){ /* col (*aj)*bs+k */ 96a7a3a9ebSHong Zhang *(bj + rowstart[itmp]) = (*aj)*bs+k; 97a7a3a9ebSHong Zhang *(bv + rowstart[itmp]) = *(av+k*bs+j); 98a7a3a9ebSHong Zhang rowstart[itmp]++; 99a7a3a9ebSHong Zhang } 100a7a3a9ebSHong Zhang } 101a7a3a9ebSHong Zhang aj++; av += bs2; 1024e5e7fe4SHong Zhang } 1034e5e7fe4SHong Zhang } 1044e5e7fe4SHong Zhang ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1054e5e7fe4SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1064e5e7fe4SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1074e5e7fe4SHong Zhang 108ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 109*c3d102feSKris Buschelman ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 110*c3d102feSKris Buschelman } else { 1114e5e7fe4SHong Zhang *newmat = B; 112*c3d102feSKris Buschelman } 1134e5e7fe4SHong Zhang PetscFunctionReturn(0); 1144e5e7fe4SHong Zhang } 1154e5e7fe4SHong Zhang #undef __FUNCT__ 11659557b74SHong Zhang #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ" 117ceb03754SKris Buschelman PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) { 118676c34cdSKris Buschelman Mat B; 11959557b74SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 120861ba921SHong Zhang Mat_SeqSBAIJ *b; 121dfbe8321SBarry Smith PetscErrorCode ierr; 12279416396SBarry Smith PetscInt *ai=a->i,*aj,m=A->M,n=A->N,i,j,*bi,*bj,*rowlengths; 123861ba921SHong Zhang PetscScalar *av,*bv; 12459557b74SHong Zhang 12559557b74SHong Zhang PetscFunctionBegin; 1262d9a3abdSHong Zhang if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 12718f87311SHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ 12859557b74SHong Zhang 12913f74950SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 13059557b74SHong Zhang for (i=0; i<m; i++) { 13159557b74SHong Zhang rowlengths[i] = ai[i+1] - a->diag[i]; 13259557b74SHong Zhang } 133f204ca49SKris Buschelman ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 134f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 135f204ca49SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr); 13659557b74SHong Zhang 137676c34cdSKris Buschelman ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 138676c34cdSKris Buschelman ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 139676c34cdSKris Buschelman ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 14059557b74SHong Zhang 141676c34cdSKris Buschelman b = (Mat_SeqSBAIJ*)(B->data); 142861ba921SHong Zhang bi = b->i; 143861ba921SHong Zhang bj = b->j; 144861ba921SHong Zhang bv = b->a; 145861ba921SHong Zhang 146861ba921SHong Zhang bi[0] = 0; 14759557b74SHong Zhang for (i=0; i<m; i++) { 14859557b74SHong Zhang aj = a->j + a->diag[i]; 14959557b74SHong Zhang av = a->a + a->diag[i]; 150861ba921SHong Zhang for (j=0; j<rowlengths[i]; j++){ 151861ba921SHong Zhang *bj = *aj; bj++; aj++; 152861ba921SHong Zhang *bv = *av; bv++; av++; 153861ba921SHong Zhang } 154861ba921SHong Zhang bi[i+1] = bi[i] + rowlengths[i]; 155861ba921SHong Zhang b->ilen[i] = rowlengths[i]; 15659557b74SHong Zhang } 15759557b74SHong Zhang 15859557b74SHong Zhang ierr = PetscFree(rowlengths);CHKERRQ(ierr); 159676c34cdSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160676c34cdSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 161676c34cdSKris Buschelman 162ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 163*c3d102feSKris Buschelman ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 164*c3d102feSKris Buschelman } else { 165676c34cdSKris Buschelman *newmat = B; 166*c3d102feSKris Buschelman } 16759557b74SHong Zhang PetscFunctionReturn(0); 16859557b74SHong Zhang } 16959557b74SHong Zhang EXTERN_C_END 17059557b74SHong Zhang 171a0e1a404SHong Zhang EXTERN_C_BEGIN 172a0e1a404SHong Zhang #undef __FUNCT__ 173a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ" 174ceb03754SKris Buschelman PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) 175a0e1a404SHong Zhang { 176a0e1a404SHong Zhang Mat B; 177a0e1a404SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 178a0e1a404SHong Zhang Mat_SeqBAIJ *b; 179dfbe8321SBarry Smith PetscErrorCode ierr; 18013f74950SBarry Smith PetscInt *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp; 18113f74950SBarry Smith PetscInt bs=A->bs,bs2=bs*bs,mbs=m/bs; 182a0e1a404SHong Zhang PetscScalar *av,*bv; 183a0e1a404SHong Zhang 184a0e1a404SHong Zhang PetscFunctionBegin; 185a0e1a404SHong Zhang /* compute browlengths of newmat */ 18613f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); 187a0e1a404SHong Zhang browstart = browlengths + mbs; 188a0e1a404SHong Zhang for (i=0; i<mbs; i++) browlengths[i] = 0; 189a0e1a404SHong Zhang aj = a->j; 190a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 191a0e1a404SHong Zhang nz = ai[i+1] - ai[i]; 192a0e1a404SHong Zhang aj++; /* skip diagonal */ 193a0e1a404SHong Zhang for (k=1; k<nz; k++) { /* no. of lower triangular blocks */ 194a0e1a404SHong Zhang browlengths[*aj]++; aj++; 195a0e1a404SHong Zhang } 196a0e1a404SHong Zhang browlengths[i] += nz; /* no. of upper triangular blocks */ 197a0e1a404SHong Zhang } 198a0e1a404SHong Zhang 199f204ca49SKris Buschelman ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 200f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 201f204ca49SKris Buschelman ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 202a0e1a404SHong Zhang ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 203a0e1a404SHong Zhang ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 204a0e1a404SHong Zhang ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 205a0e1a404SHong Zhang 206a0e1a404SHong Zhang b = (Mat_SeqBAIJ*)(B->data); 207a0e1a404SHong Zhang bi = b->i; 208a0e1a404SHong Zhang bj = b->j; 209a0e1a404SHong Zhang bv = b->a; 210a0e1a404SHong Zhang 211a0e1a404SHong Zhang /* set b->i */ 212a0e1a404SHong Zhang bi[0] = 0; 213a0e1a404SHong Zhang for (i=0; i<mbs; i++){ 214a0e1a404SHong Zhang b->ilen[i] = browlengths[i]; 215a0e1a404SHong Zhang bi[i+1] = bi[i] + browlengths[i]; 216a0e1a404SHong Zhang browstart[i] = bi[i]; 217a0e1a404SHong Zhang } 21877431f27SBarry 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); 219a0e1a404SHong Zhang 220a0e1a404SHong Zhang /* set b->j and b->a */ 221a0e1a404SHong Zhang aj = a->j; av = a->a; 222a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 223a0e1a404SHong Zhang /* diagonal block */ 224a0e1a404SHong Zhang *(bj + browstart[i]) = *aj; aj++; 225a0e1a404SHong Zhang itmp = bs2*browstart[i]; 226a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 227a0e1a404SHong Zhang *(bv + itmp + k) = *av; av++; 228a0e1a404SHong Zhang } 229a0e1a404SHong Zhang browstart[i]++; 230a0e1a404SHong Zhang 231a0e1a404SHong Zhang nz = ai[i+1] - ai[i] -1; 232a0e1a404SHong Zhang while (nz--){ 233a0e1a404SHong Zhang /* lower triangular blocks */ 234a0e1a404SHong Zhang *(bj + browstart[*aj]) = i; 235a0e1a404SHong Zhang itmp = bs2*browstart[*aj]; 236a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 237a0e1a404SHong Zhang *(bv + itmp + k) = *(av + k); 238a0e1a404SHong Zhang } 239a0e1a404SHong Zhang browstart[*aj]++; 240a0e1a404SHong Zhang 241a0e1a404SHong Zhang /* upper triangular blocks */ 242a0e1a404SHong Zhang *(bj + browstart[i]) = *aj; aj++; 243a0e1a404SHong Zhang itmp = bs2*browstart[i]; 244a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 245a0e1a404SHong Zhang *(bv + itmp + k) = *av; av++; 246a0e1a404SHong Zhang } 247a0e1a404SHong Zhang browstart[i]++; 248a0e1a404SHong Zhang } 249a0e1a404SHong Zhang } 250a0e1a404SHong Zhang ierr = PetscFree(browlengths);CHKERRQ(ierr); 251a0e1a404SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252a0e1a404SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 253a0e1a404SHong Zhang 254ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 255*c3d102feSKris Buschelman ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 256*c3d102feSKris Buschelman } else { 257a0e1a404SHong Zhang *newmat = B; 258*c3d102feSKris Buschelman } 259a0e1a404SHong Zhang PetscFunctionReturn(0); 260a0e1a404SHong Zhang } 261a0e1a404SHong Zhang #undef __FUNCT__ 262a0e1a404SHong Zhang #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ" 263ceb03754SKris Buschelman PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) 264a0e1a404SHong Zhang { 265a0e1a404SHong Zhang Mat B; 266a0e1a404SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 267a0e1a404SHong Zhang Mat_SeqSBAIJ *b; 268dfbe8321SBarry Smith PetscErrorCode ierr; 26913f74950SBarry Smith PetscInt *ai=a->i,*aj,m=A->m,n=A->n,i,j,k,*bi,*bj,*browlengths; 27013f74950SBarry Smith PetscInt bs=A->bs,bs2=bs*bs,mbs=m/bs; 271a0e1a404SHong Zhang PetscScalar *av,*bv; 272a0e1a404SHong Zhang 273a0e1a404SHong Zhang PetscFunctionBegin; 274a0e1a404SHong Zhang if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 275595208bcSHong Zhang ierr = MatMissingDiagonal_SeqBAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ 276a0e1a404SHong Zhang 27713f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); 278a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 279a0e1a404SHong Zhang browlengths[i] = ai[i+1] - a->diag[i]; 280a0e1a404SHong Zhang } 281a0e1a404SHong Zhang 282f204ca49SKris Buschelman ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 283f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 284f204ca49SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 285a0e1a404SHong Zhang ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 286a0e1a404SHong Zhang ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 287a0e1a404SHong Zhang ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 288a0e1a404SHong Zhang 289a0e1a404SHong Zhang b = (Mat_SeqSBAIJ*)(B->data); 290a0e1a404SHong Zhang bi = b->i; 291a0e1a404SHong Zhang bj = b->j; 292a0e1a404SHong Zhang bv = b->a; 293a0e1a404SHong Zhang 294a0e1a404SHong Zhang bi[0] = 0; 295a0e1a404SHong Zhang for (i=0; i<mbs; i++) { 296a0e1a404SHong Zhang aj = a->j + a->diag[i]; 297a0e1a404SHong Zhang av = a->a + (a->diag[i])*bs2; 298a0e1a404SHong Zhang for (j=0; j<browlengths[i]; j++){ 299a0e1a404SHong Zhang *bj = *aj; bj++; aj++; 300a0e1a404SHong Zhang for (k=0; k<bs2; k++){ 301a0e1a404SHong Zhang *bv = *av; bv++; av++; 302a0e1a404SHong Zhang } 303a0e1a404SHong Zhang } 304a0e1a404SHong Zhang bi[i+1] = bi[i] + browlengths[i]; 305a0e1a404SHong Zhang b->ilen[i] = browlengths[i]; 306a0e1a404SHong Zhang } 307a0e1a404SHong Zhang ierr = PetscFree(browlengths);CHKERRQ(ierr); 308a0e1a404SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 309a0e1a404SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 310a0e1a404SHong Zhang 311ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 312*c3d102feSKris Buschelman ierr = MatHeaderCopy(A);CHKERRQ(ierr); 313*c3d102feSKris Buschelman } else { 314a0e1a404SHong Zhang *newmat = B; 315*c3d102feSKris Buschelman } 316a0e1a404SHong Zhang 317a0e1a404SHong Zhang PetscFunctionReturn(0); 318a0e1a404SHong Zhang } 319a0e1a404SHong Zhang EXTERN_C_END 320