1 2 #include "src/mat/matimpl.h" 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "MatConvert_Basic" 6 /* 7 MatConvert_Basic - Converts from any input format to another format. For 8 parallel formats, the new matrix distribution is determined by PETSc. 9 10 Does not do preallocation so in general will be slow 11 */ 12 PetscErrorCode MatConvert_Basic(Mat mat,const MatType newtype,MatReuse reuse,Mat *newmat) 13 { 14 Mat M; 15 const PetscScalar *vwork; 16 PetscErrorCode ierr; 17 PetscInt i,nz,m,n,rstart,rend,lm,ln; 18 const PetscInt *cwork; 19 20 PetscFunctionBegin; 21 ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 22 ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 23 24 if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */ 25 26 ierr = MatCreate(mat->comm,lm,ln,m,n,&M);CHKERRQ(ierr); 27 ierr = MatSetType(M,newtype);CHKERRQ(ierr); 28 29 ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 30 for (i=rstart; i<rend; i++) { 31 ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 32 ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 33 ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 34 } 35 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 36 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 37 38 if (reuse == MAT_REUSE_MATRIX) { 39 ierr = MatHeaderCopy(mat,M);CHKERRQ(ierr); 40 } else { 41 *newmat = M; 42 } 43 PetscFunctionReturn(0); 44 } 45