xref: /petsc/src/mat/utils/convert.c (revision 3dab1f838b7396b990b0a386c0f17a65376f0595)
1af0996ceSBarry Smith #include <petsc/private/matimpl.h>
28b6375c0SLois Curfman McInnes 
38b6375c0SLois Curfman McInnes /*
48af18dd8SStefano Zampini   MatConvert_Basic - Converts from any input format to another format.
5273d9f13SBarry Smith   Does not do preallocation so in general will be slow
68b6375c0SLois Curfman McInnes  */
MatConvert_Basic(Mat mat,MatType newtype,MatReuse reuse,Mat * newmat)7239f794eSPierre Jolivet PetscErrorCode MatConvert_Basic(Mat mat, MatType newtype, MatReuse reuse, Mat *newmat)
8d71ae5a4SJacob Faibussowitsch {
9676c34cdSKris Buschelman   Mat                M;
10b3cc6726SBarry Smith   const PetscScalar *vwork;
112c0366f6SLisandro Dalcin   PetscInt           i, rstart, rend, nz;
12c1ac3661SBarry Smith   const PetscInt    *cwork;
1398e916f6SBarry Smith   PetscBool          isSBAIJ;
1425cdf11fSBarry Smith 
153a40ed3dSBarry Smith   PetscFunctionBegin;
168af18dd8SStefano Zampini   if (!mat->ops->getrow) { /* missing get row, use matvecs */
179566063dSJacob Faibussowitsch     PetscCall(MatConvert_Shell(mat, newtype, reuse, newmat));
183ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
198af18dd8SStefano Zampini   }
209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSBAIJ, &isSBAIJ));
2148a46eb9SPierre Jolivet   if (!isSBAIJ) PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISBAIJ, &isSBAIJ));
2228b400f6SJacob Faibussowitsch   PetscCheck(!isSBAIJ, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix");
2398e916f6SBarry Smith 
242c0366f6SLisandro Dalcin   if (reuse == MAT_REUSE_MATRIX) {
252c0366f6SLisandro Dalcin     M = *newmat;
262c0366f6SLisandro Dalcin   } else {
272c0366f6SLisandro Dalcin     PetscInt m, n, lm, ln;
289566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mat, &m, &n));
299566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mat, &lm, &ln));
309566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &M));
319566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(M, lm, ln, m, n));
329566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(M, mat, mat));
339566063dSJacob Faibussowitsch     PetscCall(MatSetType(M, newtype));
349566063dSJacob Faibussowitsch     PetscCall(MatSetUp(M));
352c0366f6SLisandro Dalcin 
369566063dSJacob Faibussowitsch     PetscCall(MatSetOption(M, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
379566063dSJacob Faibussowitsch     PetscCall(MatSetOption(M, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
38*7bccdc57SPierre Jolivet     PetscCall(MatSetOption(M, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
399566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)M, MATSEQSBAIJ, &isSBAIJ));
4048a46eb9SPierre Jolivet     if (!isSBAIJ) PetscCall(PetscObjectTypeCompare((PetscObject)M, MATMPISBAIJ, &isSBAIJ));
411baa6e33SBarry Smith     if (isSBAIJ) PetscCall(MatSetOption(M, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
422c0366f6SLisandro Dalcin   }
431a0cae60SJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(mat, &rstart, &rend));
458b6375c0SLois Curfman McInnes   for (i = rstart; i < rend; i++) {
469566063dSJacob Faibussowitsch     PetscCall(MatGetRow(mat, i, &nz, &cwork, &vwork));
479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(M, 1, &i, nz, cwork, vwork, INSERT_VALUES));
489566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(mat, i, &nz, &cwork, &vwork));
498b6375c0SLois Curfman McInnes   }
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY));
519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY));
52*7bccdc57SPierre Jolivet   if (reuse != MAT_REUSE_MATRIX) PetscCall(MatSetOption(M, MAT_NO_OFF_PROC_ENTRIES, PETSC_FALSE));
53676c34cdSKris Buschelman 
54511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
559566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(mat, &M));
56c3d102feSKris Buschelman   } else {
5718f87311SHong Zhang     *newmat = M;
58c3d102feSKris Buschelman   }
593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
608b6375c0SLois Curfman McInnes }
61