xref: /petsc/src/mat/utils/convert.c (revision c3d102fe2dfc19c1832704c537829f04ceb136a2)
156fe5c5cSLois Curfman McInnes 
294a9d846SBarry Smith #include "src/mat/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  */
12ceb03754SKris Buschelman 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;
1925cdf11fSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
218b6375c0SLois Curfman McInnes   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
22435da068SBarry Smith   ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr);
23273d9f13SBarry Smith 
24435da068SBarry Smith   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
25435da068SBarry Smith 
2618f87311SHong Zhang   ierr = MatCreate(mat->comm,lm,ln,m,n,&M);CHKERRQ(ierr);
27676c34cdSKris Buschelman   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
28273d9f13SBarry Smith 
29f1af5d2fSBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
308b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
318b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
32676c34cdSKris Buschelman     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
338b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
348b6375c0SLois Curfman McInnes   }
35676c34cdSKris Buschelman   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36676c34cdSKris Buschelman   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37676c34cdSKris Buschelman 
38ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
39*c3d102feSKris Buschelman     ierr = MatHeaderCopy(mat,M);CHKERRQ(ierr);
40*c3d102feSKris Buschelman   } else {
4118f87311SHong Zhang     *newmat = M;
42*c3d102feSKris Buschelman   }
433a40ed3dSBarry Smith   PetscFunctionReturn(0);
448b6375c0SLois Curfman McInnes }
45