156fe5c5cSLois Curfman McInnes 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> 38b6375c0SLois Curfman McInnes 48b6375c0SLois Curfman McInnes /* 58b6375c0SLois Curfman McInnes MatConvert_Basic - Converts from any input format to another format. For 68b6375c0SLois Curfman McInnes parallel formats, the new matrix distribution is determined by PETSc. 7273d9f13SBarry Smith 8273d9f13SBarry Smith Does not do preallocation so in general will be slow 98b6375c0SLois Curfman McInnes */ 10cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat mat, MatType newtype,MatReuse reuse,Mat *newmat) 11b3cc6726SBarry Smith { 12676c34cdSKris Buschelman Mat M; 13b3cc6726SBarry Smith const PetscScalar *vwork; 14dfbe8321SBarry Smith PetscErrorCode ierr; 15*2c0366f6SLisandro Dalcin PetscInt i,rstart,rend,nz; 16c1ac3661SBarry Smith const PetscInt *cwork; 1798e916f6SBarry Smith PetscBool isSBAIJ; 1825cdf11fSBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 2098e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 2198e916f6SBarry Smith if (!isSBAIJ) { 2298e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 2398e916f6SBarry Smith } 2498e916f6SBarry Smith if (isSBAIJ) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix"); 2598e916f6SBarry Smith 26*2c0366f6SLisandro Dalcin if (reuse == MAT_REUSE_MATRIX) { 27*2c0366f6SLisandro Dalcin M = *newmat; 28*2c0366f6SLisandro Dalcin } else { 29*2c0366f6SLisandro Dalcin PetscInt m,n,lm,ln; 308b6375c0SLois Curfman McInnes ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 31435da068SBarry Smith ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 32435da068SBarry Smith if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */ 33ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)mat),&M);CHKERRQ(ierr); 34f69a0ea3SMatthew Knepley ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 3533d57670SJed Brown ierr = MatSetBlockSizesFromMats(M,mat,mat);CHKERRQ(ierr); 36676c34cdSKris Buschelman ierr = MatSetType(M,newtype);CHKERRQ(ierr); 37d3bf5c86SHong Zhang ierr = MatSeqDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 38d3bf5c86SHong Zhang ierr = MatMPIDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 3998e916f6SBarry Smith ierr = MatSetUp(M);CHKERRQ(ierr); 40*2c0366f6SLisandro Dalcin 4198e916f6SBarry Smith ierr = MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 4298e916f6SBarry Smith ierr = MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 4398e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 4498e916f6SBarry Smith if (!isSBAIJ) { 4598e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 461a0cae60SJed Brown } 4798e916f6SBarry Smith if (isSBAIJ) { 4898e916f6SBarry Smith ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 49d3bf5c86SHong Zhang } 50*2c0366f6SLisandro Dalcin } 511a0cae60SJed Brown 5298e916f6SBarry Smith ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 538b6375c0SLois Curfman McInnes for (i=rstart; i<rend; i++) { 548b6375c0SLois Curfman McInnes ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 55676c34cdSKris Buschelman ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 568b6375c0SLois Curfman McInnes ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 578b6375c0SLois Curfman McInnes } 58676c34cdSKris Buschelman ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59676c34cdSKris Buschelman ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60676c34cdSKris Buschelman 61511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 6228be2f97SBarry Smith ierr = MatHeaderReplace(mat,&M);CHKERRQ(ierr); 63c3d102feSKris Buschelman } else { 6418f87311SHong Zhang *newmat = M; 65c3d102feSKris Buschelman } 663a40ed3dSBarry Smith PetscFunctionReturn(0); 678b6375c0SLois Curfman McInnes } 68