xref: /petsc/src/mat/utils/convert.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
156fe5c5cSLois Curfman McInnes 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>
38b6375c0SLois Curfman McInnes 
48b6375c0SLois Curfman McInnes /*
58af18dd8SStefano Zampini   MatConvert_Basic - Converts from any input format to another format.
6273d9f13SBarry Smith   Does not do preallocation so in general will be slow
78b6375c0SLois Curfman McInnes  */
8cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat mat,MatType newtype,MatReuse reuse,Mat *newmat)
9b3cc6726SBarry Smith {
10676c34cdSKris Buschelman   Mat               M;
11b3cc6726SBarry Smith   const PetscScalar *vwork;
122c0366f6SLisandro Dalcin   PetscInt          i,rstart,rend,nz;
13c1ac3661SBarry Smith   const PetscInt    *cwork;
1498e916f6SBarry Smith   PetscBool         isSBAIJ;
1525cdf11fSBarry Smith 
163a40ed3dSBarry Smith   PetscFunctionBegin;
178af18dd8SStefano Zampini   if (!mat->ops->getrow) { /* missing get row, use matvecs */
189566063dSJacob Faibussowitsch     PetscCall(MatConvert_Shell(mat,newtype,reuse,newmat));
198af18dd8SStefano Zampini     PetscFunctionReturn(0);
208af18dd8SStefano Zampini   }
219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&isSBAIJ));
2298e916f6SBarry Smith   if (!isSBAIJ) {
239566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&isSBAIJ));
2498e916f6SBarry Smith   }
2528b400f6SJacob Faibussowitsch   PetscCheck(!isSBAIJ,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix");
2698e916f6SBarry Smith 
272c0366f6SLisandro Dalcin   if (reuse == MAT_REUSE_MATRIX) {
282c0366f6SLisandro Dalcin     M = *newmat;
292c0366f6SLisandro Dalcin   } else {
302c0366f6SLisandro Dalcin     PetscInt m,n,lm,ln;
319566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mat,&m,&n));
329566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mat,&lm,&ln));
339566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat),&M));
349566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(M,lm,ln,m,n));
359566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(M,mat,mat));
369566063dSJacob Faibussowitsch     PetscCall(MatSetType(M,newtype));
379566063dSJacob Faibussowitsch     PetscCall(MatSetUp(M));
382c0366f6SLisandro Dalcin 
399566063dSJacob Faibussowitsch     PetscCall(MatSetOption(M,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
409566063dSJacob Faibussowitsch     PetscCall(MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
419566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isSBAIJ));
4298e916f6SBarry Smith     if (!isSBAIJ) {
439566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&isSBAIJ));
441a0cae60SJed Brown     }
45*1baa6e33SBarry Smith     if (isSBAIJ) PetscCall(MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
462c0366f6SLisandro Dalcin   }
471a0cae60SJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(mat,&rstart,&rend));
498b6375c0SLois Curfman McInnes   for (i=rstart; i<rend; i++) {
509566063dSJacob Faibussowitsch     PetscCall(MatGetRow(mat,i,&nz,&cwork,&vwork));
519566063dSJacob Faibussowitsch     PetscCall(MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES));
529566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(mat,i,&nz,&cwork,&vwork));
538b6375c0SLois Curfman McInnes   }
549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
56676c34cdSKris Buschelman 
57511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
589566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(mat,&M));
59c3d102feSKris Buschelman   } else {
6018f87311SHong Zhang     *newmat = M;
61c3d102feSKris Buschelman   }
623a40ed3dSBarry Smith   PetscFunctionReturn(0);
638b6375c0SLois Curfman McInnes }
64