1 2 #include <petsc/private/matimpl.h> 3 4 /* 5 MatConvert_Basic - Converts from any input format to another format. 6 Does not do preallocation so in general will be slow 7 */ 8 PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat mat,MatType newtype,MatReuse reuse,Mat *newmat) 9 { 10 Mat M; 11 const PetscScalar *vwork; 12 PetscErrorCode ierr; 13 PetscInt i,rstart,rend,nz; 14 const PetscInt *cwork; 15 PetscBool isSBAIJ; 16 17 PetscFunctionBegin; 18 if (!mat->ops->getrow) { /* missing get row, use matvecs */ 19 ierr = MatConvert_Shell(mat,newtype,reuse,newmat);CHKERRQ(ierr); 20 PetscFunctionReturn(0); 21 } 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 if (reuse == MAT_REUSE_MATRIX) { 29 M = *newmat; 30 } else { 31 PetscInt m,n,lm,ln; 32 ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 33 ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr); 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 ierr = MatSeqDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 39 ierr = MatMPIDenseSetPreallocation(M,NULL);CHKERRQ(ierr); 40 ierr = MatSetUp(M);CHKERRQ(ierr); 41 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 ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr); 45 if (!isSBAIJ) { 46 ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr); 47 } 48 if (isSBAIJ) { 49 ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 50 } 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_INPLACE_MATRIX) { 63 ierr = MatHeaderReplace(mat,&M);CHKERRQ(ierr); 64 } else { 65 *newmat = M; 66 } 67 PetscFunctionReturn(0); 68 } 69