Actual source code: convert.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2:  #include <petsc/private/matimpl.h>

  4: /*
  5:   MatConvert_Basic - Converts from any input format to another format. For
  6:   parallel formats, the new matrix distribution is determined by PETSc.

  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;

 20:   PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ);
 21:   if (!isSBAIJ) {
 22:     PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ);
 23:   }
 24:   if (isSBAIJ) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix");

 26:   if (reuse == MAT_REUSE_MATRIX) {
 27:     M = *newmat;
 28:   } else {
 29:     PetscInt m,n,lm,ln;
 30:     MatGetSize(mat,&m,&n);
 31:     MatGetLocalSize(mat,&lm,&ln);
 32:     if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
 33:     MatCreate(PetscObjectComm((PetscObject)mat),&M);
 34:     MatSetSizes(M,lm,ln,m,n);
 35:     MatSetBlockSizesFromMats(M,mat,mat);
 36:     MatSetType(M,newtype);
 37:     MatSeqDenseSetPreallocation(M,NULL);
 38:     MatMPIDenseSetPreallocation(M,NULL);
 39:     MatSetUp(M);

 41:     MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
 42:     MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
 43:     PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);
 44:     if (!isSBAIJ) {
 45:       PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);
 46:     }
 47:     if (isSBAIJ) {
 48:       MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
 49:     }
 50:   }

 52:   MatGetOwnershipRange(mat,&rstart,&rend);
 53:   for (i=rstart; i<rend; i++) {
 54:     MatGetRow(mat,i,&nz,&cwork,&vwork);
 55:     MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);
 56:     MatRestoreRow(mat,i,&nz,&cwork,&vwork);
 57:   }
 58:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 61:   if (reuse == MAT_INPLACE_MATRIX) {
 62:     MatHeaderReplace(mat,&M);
 63:   } else {
 64:     *newmat = M;
 65:   }
 66:   return(0);
 67: }