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