1 2 #include <petsc/private/matimpl.h> 3 4 /* 5 MatConvert_Basic - Converts from any input format to another format. For 6 parallel formats, the new matrix distribution is determined by PETSc. 7 8 Does not do preallocation so in general will be slow 9 */ 10 PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat mat, MatType newtype,MatReuse reuse,Mat *newmat) 11 { 12 Mat M; 13 const PetscScalar *vwork; 14 PetscErrorCode ierr; 15 PetscInt i,rstart,rend,nz; 16 const PetscInt *cwork; 17 PetscBool isSBAIJ; 18 19 PetscFunctionBegin; 20 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 21 if (!isSBAIJ) { 22 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 23 } 24 if (isSBAIJ) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix"); 25 26 if (reuse == MAT_REUSE_MATRIX) { 27 M = *newmat; 28 } else { 29 PetscInt m,n,lm,ln; 30 ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 31 ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 32 if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */ 33 ierr = MatCreate(PetscObjectComm((PetscObject)mat),&M);CHKERRQ(ierr); 34 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 35 ierr = MatSetBlockSizesFromMats(M,mat,mat);CHKERRQ(ierr); 36 ierr = MatSetType(M,newtype);CHKERRQ(ierr); 37 ierr = MatSeqDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 38 ierr = MatMPIDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 39 ierr = MatSetUp(M);CHKERRQ(ierr); 40 41 ierr = MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 42 ierr = MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 43 ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 44 if (!isSBAIJ) { 45 ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 46 } 47 if (isSBAIJ) { 48 ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 49 } 50 } 51 52 ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 53 for (i=rstart; i<rend; i++) { 54 ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 55 ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 56 ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 57 } 58 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60 61 if (reuse == MAT_INPLACE_MATRIX) { 62 ierr = MatHeaderReplace(mat,&M);CHKERRQ(ierr); 63 } else { 64 *newmat = M; 65 } 66 PetscFunctionReturn(0); 67 } 68