xref: /petsc/src/mat/utils/convert.c (revision 251f4c6745e69799cb0300fd80598c58f96eb452)
156fe5c5cSLois Curfman McInnes 
2b45d2f2cSJed Brown #include <petsc-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);
29535b19f3SBarry Smith   ierr = MatSetBlockSize(M,mat->rmap->bs);CHKERRQ(ierr);
30676c34cdSKris Buschelman   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
31*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
32*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
334994cf47SJed Brown   ierr = MatSetUp(M);CHKERRQ(ierr);
34d1bf975fSJed Brown   if (isseqsbaij || ismpisbaij) {ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);}
35273d9f13SBarry Smith 
36f1af5d2fSBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
378b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
388b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
39676c34cdSKris Buschelman     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
408b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
418b6375c0SLois Curfman McInnes   }
42676c34cdSKris Buschelman   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
43676c34cdSKris Buschelman   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
44676c34cdSKris Buschelman 
45ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
468ab5b326SKris Buschelman     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
47c3d102feSKris Buschelman   } else {
4818f87311SHong Zhang     *newmat = M;
49c3d102feSKris Buschelman   }
503a40ed3dSBarry Smith   PetscFunctionReturn(0);
518b6375c0SLois Curfman McInnes }
52