Actual source code: convert.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: #include <petsc/private/matimpl.h>

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

 10:   Does not do preallocation so in general will be slow
 11:  */
 12: PETSC_INTERN 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;

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


 29:   MatGetSize(mat,&m,&n);
 30:   MatGetLocalSize(mat,&lm,&ln);

 32:   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */

 34:   MatCreate(PetscObjectComm((PetscObject)mat),&M);
 35:   MatSetSizes(M,lm,ln,m,n);
 36:   MatSetBlockSizesFromMats(M,mat,mat);
 37:   MatSetType(M,newtype);

 39:   MatSeqDenseSetPreallocation(M,NULL);
 40:   MatMPIDenseSetPreallocation(M,NULL);
 41:   MatSetUp(M);
 42:   MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
 43:   MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

 45:     PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);
 46:   if (!isSBAIJ) {
 47:     PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);
 48:   }
 49:   if (isSBAIJ) {
 50:     MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
 51:   }

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

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