156fe5c5cSLois Curfman McInnes 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> 38b6375c0SLois Curfman McInnes 48b6375c0SLois Curfman McInnes /* 58af18dd8SStefano Zampini MatConvert_Basic - Converts from any input format to another format. 6273d9f13SBarry Smith Does not do preallocation so in general will be slow 78b6375c0SLois Curfman McInnes */ 8cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat mat,MatType newtype,MatReuse reuse,Mat *newmat) 9b3cc6726SBarry Smith { 10676c34cdSKris Buschelman Mat M; 11b3cc6726SBarry Smith const PetscScalar *vwork; 12dfbe8321SBarry Smith PetscErrorCode ierr; 132c0366f6SLisandro Dalcin PetscInt i,rstart,rend,nz; 14c1ac3661SBarry Smith const PetscInt *cwork; 1598e916f6SBarry Smith PetscBool isSBAIJ; 1625cdf11fSBarry Smith 173a40ed3dSBarry Smith PetscFunctionBegin; 188af18dd8SStefano Zampini if (!mat->ops->getrow) { /* missing get row, use matvecs */ 198af18dd8SStefano Zampini ierr = MatConvert_Shell(mat,newtype,reuse,newmat);CHKERRQ(ierr); 208af18dd8SStefano Zampini PetscFunctionReturn(0); 218af18dd8SStefano Zampini } 2298e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 2398e916f6SBarry Smith if (!isSBAIJ) { 2498e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 2598e916f6SBarry Smith } 26*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(isSBAIJ,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix"); 2798e916f6SBarry Smith 282c0366f6SLisandro Dalcin if (reuse == MAT_REUSE_MATRIX) { 292c0366f6SLisandro Dalcin M = *newmat; 302c0366f6SLisandro Dalcin } else { 312c0366f6SLisandro Dalcin PetscInt m,n,lm,ln; 328b6375c0SLois Curfman McInnes ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 33435da068SBarry Smith ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 34ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)mat),&M);CHKERRQ(ierr); 35f69a0ea3SMatthew Knepley ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 3633d57670SJed Brown ierr = MatSetBlockSizesFromMats(M,mat,mat);CHKERRQ(ierr); 37676c34cdSKris Buschelman ierr = MatSetType(M,newtype);CHKERRQ(ierr); 38d3bf5c86SHong Zhang ierr = MatSeqDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 39d3bf5c86SHong Zhang ierr = MatMPIDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 4098e916f6SBarry Smith ierr = MatSetUp(M);CHKERRQ(ierr); 412c0366f6SLisandro Dalcin 4298e916f6SBarry Smith ierr = MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 4398e916f6SBarry Smith ierr = MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 4498e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 4598e916f6SBarry Smith if (!isSBAIJ) { 4698e916f6SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 471a0cae60SJed Brown } 4898e916f6SBarry Smith if (isSBAIJ) { 4998e916f6SBarry Smith ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 50d3bf5c86SHong Zhang } 512c0366f6SLisandro Dalcin } 521a0cae60SJed Brown 5398e916f6SBarry Smith ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 548b6375c0SLois Curfman McInnes for (i=rstart; i<rend; i++) { 558b6375c0SLois Curfman McInnes ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 56676c34cdSKris Buschelman ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 578b6375c0SLois Curfman McInnes ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 588b6375c0SLois Curfman McInnes } 59676c34cdSKris Buschelman ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60676c34cdSKris Buschelman ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61676c34cdSKris Buschelman 62511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 6328be2f97SBarry Smith ierr = MatHeaderReplace(mat,&M);CHKERRQ(ierr); 64c3d102feSKris Buschelman } else { 6518f87311SHong Zhang *newmat = M; 66c3d102feSKris Buschelman } 673a40ed3dSBarry Smith PetscFunctionReturn(0); 688b6375c0SLois Curfman McInnes } 69