1*be1d678aSKris Buschelman #define PETSCMAT_DLL 256fe5c5cSLois Curfman McInnes 394a9d846SBarry Smith #include "src/mat/matimpl.h" 48b6375c0SLois Curfman McInnes 54a2ae208SSatish Balay #undef __FUNCT__ 64a2ae208SSatish Balay #define __FUNCT__ "MatConvert_Basic" 78b6375c0SLois Curfman McInnes /* 88b6375c0SLois Curfman McInnes MatConvert_Basic - Converts from any input format to another format. For 98b6375c0SLois Curfman McInnes parallel formats, the new matrix distribution is determined by PETSc. 10273d9f13SBarry Smith 11273d9f13SBarry Smith Does not do preallocation so in general will be slow 128b6375c0SLois Curfman McInnes */ 13ceb03754SKris Buschelman PetscErrorCode MatConvert_Basic(Mat mat,const MatType newtype,MatReuse reuse,Mat *newmat) 14b3cc6726SBarry Smith { 15676c34cdSKris Buschelman Mat M; 16b3cc6726SBarry Smith const PetscScalar *vwork; 17dfbe8321SBarry Smith PetscErrorCode ierr; 18c1ac3661SBarry Smith PetscInt i,nz,m,n,rstart,rend,lm,ln; 19c1ac3661SBarry Smith const PetscInt *cwork; 2025cdf11fSBarry Smith 213a40ed3dSBarry Smith PetscFunctionBegin; 228b6375c0SLois Curfman McInnes ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 23435da068SBarry Smith ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 24273d9f13SBarry Smith 25435da068SBarry Smith if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */ 26435da068SBarry Smith 2718f87311SHong Zhang ierr = MatCreate(mat->comm,lm,ln,m,n,&M);CHKERRQ(ierr); 28676c34cdSKris Buschelman ierr = MatSetType(M,newtype);CHKERRQ(ierr); 29273d9f13SBarry Smith 30f1af5d2fSBarry Smith ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 318b6375c0SLois Curfman McInnes for (i=rstart; i<rend; i++) { 328b6375c0SLois Curfman McInnes ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 33676c34cdSKris Buschelman ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 348b6375c0SLois Curfman McInnes ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 358b6375c0SLois Curfman McInnes } 36676c34cdSKris Buschelman ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 37676c34cdSKris Buschelman ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 38676c34cdSKris Buschelman 39ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 408ab5b326SKris Buschelman ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr); 41c3d102feSKris Buschelman } else { 4218f87311SHong Zhang *newmat = M; 43c3d102feSKris Buschelman } 443a40ed3dSBarry Smith PetscFunctionReturn(0); 458b6375c0SLois Curfman McInnes } 46