1 2 #include <petsc/private/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, MatType newtype,MatReuse reuse,Mat *newmat) 13 { 14 Mat M; 15 const PetscScalar *vwork; 16 PetscErrorCode ierr; 17 PetscInt nz,i,m,n,rstart,rend,lm,ln; 18 const PetscInt *cwork; 19 PetscBool isSBAIJ; 20 21 PetscFunctionBegin; 22 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 23 if (!isSBAIJ) { 24 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 25 } 26 if (isSBAIJ) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix"); 27 28 29 ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 30 ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 31 32 if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */ 33 34 ierr = MatCreate(PetscObjectComm((PetscObject)mat),&M);CHKERRQ(ierr); 35 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 36 ierr = MatSetBlockSizesFromMats(M,mat,mat);CHKERRQ(ierr); 37 ierr = MatSetType(M,newtype);CHKERRQ(ierr); 38 39 ierr = MatSeqDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 40 ierr = MatMPIDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 41 ierr = MatSetUp(M);CHKERRQ(ierr); 42 ierr = MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 43 ierr = MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 44 45 ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 46 if (!isSBAIJ) { 47 ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 48 } 49 if (isSBAIJ) { 50 ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 51 } 52 53 ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 54 for (i=rstart; i<rend; i++) { 55 ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 56 ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 57 ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 58 } 59 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61 62 if (reuse == MAT_REUSE_MATRIX) { 63 ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr); 64 } else { 65 *newmat = M; 66 } 67 PetscFunctionReturn(0); 68 } 69