xref: /petsc/src/mat/utils/convert.c (revision 4994cf479fe3cab9bbc99db3d00c8198cfe9191f)
156fe5c5cSLois Curfman McInnes 
2c6db04a5SJed Brown #include <private/matimpl.h>
38b6375c0SLois Curfman McInnes 
44a2ae208SSatish Balay #undef __FUNCT__
54a2ae208SSatish Balay #define __FUNCT__ "MatConvert_Basic"
68b6375c0SLois Curfman McInnes /*
78b6375c0SLois Curfman McInnes   MatConvert_Basic - Converts from any input format to another format. For
88b6375c0SLois Curfman McInnes   parallel formats, the new matrix distribution is determined by PETSc.
9273d9f13SBarry Smith 
10273d9f13SBarry Smith   Does not do preallocation so in general will be slow
118b6375c0SLois Curfman McInnes  */
12a313700dSBarry Smith PetscErrorCode MatConvert_Basic(Mat mat, const MatType newtype,MatReuse reuse,Mat *newmat)
13b3cc6726SBarry Smith {
14676c34cdSKris Buschelman   Mat                M;
15b3cc6726SBarry Smith   const PetscScalar  *vwork;
16dfbe8321SBarry Smith   PetscErrorCode     ierr;
17c1ac3661SBarry Smith   PetscInt           i,nz,m,n,rstart,rend,lm,ln;
18c1ac3661SBarry Smith   const PetscInt     *cwork;
19d1bf975fSJed Brown   PetscBool          isseqsbaij,ismpisbaij;
2025cdf11fSBarry Smith 
213a40ed3dSBarry Smith   PetscFunctionBegin;
228b6375c0SLois Curfman McInnes   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
23435da068SBarry Smith   ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr);
24273d9f13SBarry Smith 
25435da068SBarry Smith   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
26435da068SBarry Smith 
277adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)mat)->comm,&M);CHKERRQ(ierr);
28f69a0ea3SMatthew Knepley   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
29676c34cdSKris Buschelman   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
30d1bf975fSJed Brown   ierr = PetscTypeCompare((PetscObject)M,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
31d1bf975fSJed Brown   ierr = PetscTypeCompare((PetscObject)M,MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
32*4994cf47SJed Brown   ierr = MatSetUp(M);CHKERRQ(ierr);
33d1bf975fSJed Brown   if (isseqsbaij || ismpisbaij) {ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);}
34273d9f13SBarry Smith 
35f1af5d2fSBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
368b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
378b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
38676c34cdSKris Buschelman     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
398b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
408b6375c0SLois Curfman McInnes   }
41676c34cdSKris Buschelman   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42676c34cdSKris Buschelman   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
43676c34cdSKris Buschelman 
44ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
458ab5b326SKris Buschelman     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
46c3d102feSKris Buschelman   } else {
4718f87311SHong Zhang     *newmat = M;
48c3d102feSKris Buschelman   }
493a40ed3dSBarry Smith   PetscFunctionReturn(0);
508b6375c0SLois Curfman McInnes }
51