xref: /petsc/src/mat/utils/convert.c (revision c6db04a5321582041def2b1e244c75985478b3ef)
156fe5c5cSLois Curfman McInnes 
2*c6db04a5SJed 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;
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 
267adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)mat)->comm,&M);CHKERRQ(ierr);
27f69a0ea3SMatthew Knepley   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
28676c34cdSKris Buschelman   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
29273d9f13SBarry Smith 
30f1af5d2fSBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
318b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
328b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
33676c34cdSKris Buschelman     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
348b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
358b6375c0SLois Curfman McInnes   }
36676c34cdSKris Buschelman   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37676c34cdSKris Buschelman   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38547795f9SHong Zhang   if (mat->hermitian){
39547795f9SHong Zhang     ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
40547795f9SHong Zhang   }
41676c34cdSKris Buschelman 
42ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
438ab5b326SKris Buschelman     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
44c3d102feSKris Buschelman   } else {
4518f87311SHong Zhang     *newmat = M;
46c3d102feSKris Buschelman   }
473a40ed3dSBarry Smith   PetscFunctionReturn(0);
488b6375c0SLois Curfman McInnes }
49