xref: /petsc/src/mat/utils/convert.c (revision 98e916f64b3a11a4188c92a91329d7e81d10bedd)
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  */
1219fd82e9SBarry Smith PetscErrorCode MatConvert_Basic(Mat mat, MatType newtype,MatReuse reuse,Mat *newmat)
13b3cc6726SBarry Smith {
14676c34cdSKris Buschelman   Mat               M;
15b3cc6726SBarry Smith   const PetscScalar *vwork;
16dfbe8321SBarry Smith   PetscErrorCode    ierr;
17*98e916f6SBarry Smith   PetscInt          nz,i,m,n,rstart,rend,lm,ln;
18c1ac3661SBarry Smith   const PetscInt    *cwork;
19*98e916f6SBarry Smith   PetscBool         isSBAIJ;
2025cdf11fSBarry Smith 
213a40ed3dSBarry Smith   PetscFunctionBegin;
22*98e916f6SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr);
23*98e916f6SBarry Smith   if (!isSBAIJ) {
24*98e916f6SBarry Smith     ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr);
25*98e916f6SBarry Smith   }
26*98e916f6SBarry Smith   if (isSBAIJ) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix");
27*98e916f6SBarry Smith 
28*98e916f6SBarry Smith 
298b6375c0SLois Curfman McInnes   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
30435da068SBarry Smith   ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr);
31273d9f13SBarry Smith 
32435da068SBarry Smith   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
33435da068SBarry Smith 
34ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)mat),&M);CHKERRQ(ierr);
35f69a0ea3SMatthew Knepley   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
3633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(M,mat,mat);CHKERRQ(ierr);
37676c34cdSKris Buschelman   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
38d3bf5c86SHong Zhang 
39d3bf5c86SHong Zhang   ierr = MatSeqDenseSetPreallocation(M,NULL);CHKERRQ(ierr);
40d3bf5c86SHong Zhang   ierr = MatMPIDenseSetPreallocation(M,NULL);CHKERRQ(ierr);
41*98e916f6SBarry Smith   ierr = MatSetUp(M);CHKERRQ(ierr);
42*98e916f6SBarry Smith   ierr = MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
43*98e916f6SBarry Smith   ierr = MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
441a0cae60SJed Brown 
45*98e916f6SBarry Smith     ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ);CHKERRQ(ierr);
46*98e916f6SBarry Smith   if (!isSBAIJ) {
47*98e916f6SBarry Smith     ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ);CHKERRQ(ierr);
481a0cae60SJed Brown   }
49*98e916f6SBarry Smith   if (isSBAIJ) {
50*98e916f6SBarry Smith     ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
51d3bf5c86SHong Zhang   }
521a0cae60SJed Brown 
53*98e916f6SBarry Smith   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
548b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
558b6375c0SLois Curfman McInnes     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
56676c34cdSKris Buschelman     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
578b6375c0SLois Curfman McInnes     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
588b6375c0SLois Curfman McInnes   }
59676c34cdSKris Buschelman   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60676c34cdSKris Buschelman   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61676c34cdSKris Buschelman 
62ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
638ab5b326SKris Buschelman     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
64c3d102feSKris Buschelman   } else {
6518f87311SHong Zhang     *newmat = M;
66c3d102feSKris Buschelman   }
673a40ed3dSBarry Smith   PetscFunctionReturn(0);
688b6375c0SLois Curfman McInnes }
69