xref: /petsc/src/mat/utils/convert.c (revision 547795f9f637a72c8a9fe35eb9ea2bca0ad69cfa)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
256fe5c5cSLois Curfman McInnes 
37c4f633dSBarry Smith #include "private/matimpl.h"
48b6375c0SLois Curfman McInnes 
54a2ae208SSatish Balay #undef __FUNCT__
64a2ae208SSatish Balay #define __FUNCT__ "MatConvert_Basic"
78b6375c0SLois Curfman McInnes /*
88b6375c0SLois Curfman McInnes   MatConvert_Basic - Converts from any input format to another format. For
98b6375c0SLois Curfman McInnes   parallel formats, the new matrix distribution is determined by PETSc.
10273d9f13SBarry Smith 
11273d9f13SBarry Smith   Does not do preallocation so in general will be slow
128b6375c0SLois Curfman McInnes  */
13a313700dSBarry Smith PetscErrorCode MatConvert_Basic(Mat mat, const MatType newtype,MatReuse reuse,Mat *newmat)
14b3cc6726SBarry Smith {
15676c34cdSKris Buschelman   Mat                M;
16b3cc6726SBarry Smith   const PetscScalar  *vwork;
17dfbe8321SBarry Smith   PetscErrorCode     ierr;
18c1ac3661SBarry Smith   PetscInt           i,nz,m,n,rstart,rend,lm,ln;
19c1ac3661SBarry Smith   const PetscInt     *cwork;
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);
30273d9f13SBarry Smith 
31f1af5d2fSBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
328b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
338b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
34676c34cdSKris Buschelman     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
358b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
368b6375c0SLois Curfman McInnes   }
37676c34cdSKris Buschelman   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38676c34cdSKris Buschelman   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39*547795f9SHong Zhang   if (mat->hermitian){
40*547795f9SHong Zhang     ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
41*547795f9SHong Zhang   }
42676c34cdSKris Buschelman 
43ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
448ab5b326SKris Buschelman     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
45c3d102feSKris Buschelman   } else {
4618f87311SHong Zhang     *newmat = M;
47c3d102feSKris Buschelman   }
483a40ed3dSBarry Smith   PetscFunctionReturn(0);
498b6375c0SLois Curfman McInnes }
50