xref: /petsc/src/mat/utils/convert.c (revision 435da068ed131cbb9cf235815d431601925f48e7)
1*435da068SBarry Smith /*$Id: convert.c,v 1.72 2001/02/13 16:24:58 bsmith Exp bsmith $*/
256fe5c5cSLois Curfman McInnes 
394a9d846SBarry Smith #include "src/mat/matimpl.h"
48b6375c0SLois Curfman McInnes 
55615d1e5SSatish Balay #undef __FUNC__
6b0a32e0cSBarry Smith #define __FUNC__ "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  */
138b6375c0SLois Curfman McInnes int MatConvert_Basic(Mat mat,MatType newtype,Mat *M)
148b6375c0SLois Curfman McInnes {
158b6375c0SLois Curfman McInnes   Scalar     *vwork;
16*435da068SBarry Smith   int        ierr,i,nz,m,n,*cwork,rstart,rend,lm,ln;
1725cdf11fSBarry Smith 
183a40ed3dSBarry Smith   PetscFunctionBegin;
198b6375c0SLois Curfman McInnes   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
20*435da068SBarry Smith   ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr);
21273d9f13SBarry Smith 
22*435da068SBarry Smith   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
23*435da068SBarry Smith 
24*435da068SBarry Smith   ierr = MatCreate(mat->comm,lm,ln,m,n,M);CHKERRQ(ierr);
25273d9f13SBarry Smith   ierr = MatSetType(*M,newtype);CHKERRQ(ierr);
26273d9f13SBarry Smith 
27f1af5d2fSBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
288b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
298b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
308b6375c0SLois Curfman McInnes     ierr = MatSetValues(*M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
318b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
328b6375c0SLois Curfman McInnes   }
336d4a8577SBarry Smith   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
346d4a8577SBarry Smith   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
353a40ed3dSBarry Smith   PetscFunctionReturn(0);
368b6375c0SLois Curfman McInnes }
37