xref: /petsc/src/mat/utils/convert.c (revision 047240e14af00aad1ef65e96f6fface8924f7f7e)
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          i,nz,m,n,rstart,rend,lm,ln;
18   const PetscInt    *cwork;
19   PetscBool         isseqsbaij,ismpisbaij;
20 
21   PetscFunctionBegin;
22   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
23   ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr);
24 
25   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
26 
27   ierr = MatCreate(((PetscObject)mat)->comm,&M);CHKERRQ(ierr);
28   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
29   ierr = MatSetBlockSize(M,mat->rmap->bs);CHKERRQ(ierr);
30   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
31   ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
32   ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
33   ierr = MatSetUp(M);CHKERRQ(ierr);
34   if (isseqsbaij || ismpisbaij) {ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);}
35 
36   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
37   for (i=rstart; i<rend; i++) {
38     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
39     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
40     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
41   }
42   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
43   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
44 
45   if (reuse == MAT_REUSE_MATRIX) {
46     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
47   } else {
48     *newmat = M;
49   }
50   PetscFunctionReturn(0);
51 }
52