xref: /petsc/src/mat/utils/convert.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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  */
8d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat mat, MatType newtype, MatReuse reuse, Mat *newmat)
9d71ae5a4SJacob Faibussowitsch {
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));
19*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
208af18dd8SStefano Zampini   }
219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSBAIJ, &isSBAIJ));
2248a46eb9SPierre Jolivet   if (!isSBAIJ) PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPISBAIJ, &isSBAIJ));
2328b400f6SJacob Faibussowitsch   PetscCheck(!isSBAIJ, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot convert from SBAIJ matrix since cannot obtain entire rows of matrix");
2498e916f6SBarry Smith 
252c0366f6SLisandro Dalcin   if (reuse == MAT_REUSE_MATRIX) {
262c0366f6SLisandro Dalcin     M = *newmat;
272c0366f6SLisandro Dalcin   } else {
282c0366f6SLisandro Dalcin     PetscInt m, n, lm, ln;
299566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mat, &m, &n));
309566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mat, &lm, &ln));
319566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &M));
329566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(M, lm, ln, m, n));
339566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(M, mat, mat));
349566063dSJacob Faibussowitsch     PetscCall(MatSetType(M, newtype));
359566063dSJacob Faibussowitsch     PetscCall(MatSetUp(M));
362c0366f6SLisandro Dalcin 
379566063dSJacob Faibussowitsch     PetscCall(MatSetOption(M, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
389566063dSJacob Faibussowitsch     PetscCall(MatSetOption(M, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
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));
52676c34cdSKris Buschelman 
53511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
549566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(mat, &M));
55c3d102feSKris Buschelman   } else {
5618f87311SHong Zhang     *newmat = M;
57c3d102feSKris Buschelman   }
58*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
598b6375c0SLois Curfman McInnes }
60