xref: /petsc/src/mat/utils/axpy.c (revision b20dac60dd09ad4261f883244468bb329c680de2)
1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
26f79c3a4SBarry Smith 
3d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTransposeAXPY_Private(Mat Y, PetscScalar a, Mat X, MatStructure str, Mat T)
4d71ae5a4SJacob Faibussowitsch {
5d54c938cSPierre Jolivet   Mat         A, F;
64166f22eSPierre Jolivet   PetscScalar vshift, vscale;
75f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat, Mat *);
8d54c938cSPierre Jolivet 
9d54c938cSPierre Jolivet   PetscFunctionBegin;
104166f22eSPierre Jolivet   if (T == X) PetscCall(MatShellGetScalingShifts(T, &vshift, &vscale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
114166f22eSPierre Jolivet   else {
124166f22eSPierre Jolivet     vshift = 0.0;
134166f22eSPierre Jolivet     vscale = 1.0;
144166f22eSPierre Jolivet   }
159566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)T, "MatTransposeGetMat_C", &f));
16d54c938cSPierre Jolivet   if (f) {
179566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(T, &A));
18d54c938cSPierre Jolivet     if (T == X) {
19013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Explicitly transposing X of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
209566063dSJacob Faibussowitsch       PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &F));
21d54c938cSPierre Jolivet       A = Y;
22d54c938cSPierre Jolivet     } else {
23013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Transposing X because Y of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
249566063dSJacob Faibussowitsch       PetscCall(MatTranspose(X, MAT_INITIAL_MATRIX, &F));
25d54c938cSPierre Jolivet     }
26d54c938cSPierre Jolivet   } else {
279566063dSJacob Faibussowitsch     PetscCall(MatHermitianTransposeGetMat(T, &A));
28d54c938cSPierre Jolivet     if (T == X) {
2952cd20c4SPierre Jolivet       PetscCall(PetscInfo(NULL, "Explicitly Hermitian transposing X of type MATHERMITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
309566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &F));
31d54c938cSPierre Jolivet       A = Y;
32d54c938cSPierre Jolivet     } else {
3352cd20c4SPierre Jolivet       PetscCall(PetscInfo(NULL, "Hermitian transposing X because Y of type MATHERMITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
349566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(X, MAT_INITIAL_MATRIX, &F));
35d54c938cSPierre Jolivet     }
36d54c938cSPierre Jolivet   }
374166f22eSPierre Jolivet   PetscCall(MatAXPY(A, a * vscale, F, str));
384166f22eSPierre Jolivet   PetscCall(MatShift(A, a * vshift));
399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41d54c938cSPierre Jolivet }
42d54c938cSPierre Jolivet 
43*b20dac60SPierre Jolivet static PetscErrorCode MatAXPY_BasicWithTypeCompare(Mat Y, PetscScalar a, Mat X, MatStructure str)
44*b20dac60SPierre Jolivet {
45*b20dac60SPierre Jolivet   PetscBool flg;
46*b20dac60SPierre Jolivet 
47*b20dac60SPierre Jolivet   PetscFunctionBegin;
48*b20dac60SPierre Jolivet   PetscCall(MatIsShell(Y, &flg));
49*b20dac60SPierre Jolivet   if (flg) { /* MatShell has special support for AXPY */
50*b20dac60SPierre Jolivet     PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure);
51*b20dac60SPierre Jolivet 
52*b20dac60SPierre Jolivet     PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f));
53*b20dac60SPierre Jolivet     if (f) {
54*b20dac60SPierre Jolivet       PetscCall((*f)(Y, a, X, str));
55*b20dac60SPierre Jolivet       PetscFunctionReturn(PETSC_SUCCESS);
56*b20dac60SPierre Jolivet     }
57*b20dac60SPierre Jolivet   } else {
58*b20dac60SPierre Jolivet     /* no need to preallocate if Y is dense */
59*b20dac60SPierre Jolivet     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &flg, MATSEQDENSE, MATMPIDENSE, ""));
60*b20dac60SPierre Jolivet     if (flg) {
61*b20dac60SPierre Jolivet       PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &flg));
62*b20dac60SPierre Jolivet       if (flg) {
63*b20dac60SPierre Jolivet         PetscCall(MatAXPY_Dense_Nest(Y, a, X));
64*b20dac60SPierre Jolivet         PetscFunctionReturn(PETSC_SUCCESS);
65*b20dac60SPierre Jolivet       }
66*b20dac60SPierre Jolivet     }
67*b20dac60SPierre Jolivet     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSCALAPACK, MATELEMENTAL, ""));
68*b20dac60SPierre Jolivet     if (flg) { /* Avoid MatAXPY_Basic() due to missing MatGetRow() */
69*b20dac60SPierre Jolivet       Mat C;
70*b20dac60SPierre Jolivet 
71*b20dac60SPierre Jolivet       PetscCall(MatConvert(X, ((PetscObject)Y)->type_name, MAT_INITIAL_MATRIX, &C));
72*b20dac60SPierre Jolivet       PetscCall(MatAXPY(Y, a, C, str));
73*b20dac60SPierre Jolivet       PetscCall(MatDestroy(&C));
74*b20dac60SPierre Jolivet       PetscFunctionReturn(PETSC_SUCCESS);
75*b20dac60SPierre Jolivet     }
76*b20dac60SPierre Jolivet   }
77*b20dac60SPierre Jolivet   PetscCall(MatAXPY_Basic(Y, a, X, str));
78*b20dac60SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
79*b20dac60SPierre Jolivet }
80*b20dac60SPierre Jolivet 
8106be10caSBarry Smith /*@
8221c89e3eSBarry Smith   MatAXPY - Computes Y = a*X + Y.
836f79c3a4SBarry Smith 
84c3339decSBarry Smith   Logically Collective
85fee21e36SBarry Smith 
8698a79cdbSBarry Smith   Input Parameters:
87607cd303SBarry Smith + a   - the scalar multiplier
88607cd303SBarry Smith . X   - the first matrix
89607cd303SBarry Smith . Y   - the second matrix
9027430b45SBarry Smith - str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN`, `UNKNOWN_NONZERO_PATTERN`, or `SUBSET_NONZERO_PATTERN` (nonzeros of `X` is a subset of `Y`'s)
91edf9a46cSStefano Zampini 
922860a424SLois Curfman McInnes   Level: intermediate
932860a424SLois Curfman McInnes 
941cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAYPX()`
9506be10caSBarry Smith  @*/
96d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str)
97d71ae5a4SJacob Faibussowitsch {
98646531bbSStefano Zampini   PetscInt  M1, M2, N1, N2;
99c1ac3661SBarry Smith   PetscInt  m1, m2, n1, n2;
100*b20dac60SPierre Jolivet   PetscBool sametype, transpose;
1016f79c3a4SBarry Smith 
1023a40ed3dSBarry Smith   PetscFunctionBegin;
1030700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
104c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y, a, 2);
105646531bbSStefano Zampini   PetscValidHeaderSpecific(X, MAT_CLASSID, 3);
106646531bbSStefano Zampini   PetscCheckSameComm(Y, 1, X, 3);
1079566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M1, &N1));
1089566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M2, &N2));
1099566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m1, &n1));
1109566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m2, &n2));
1112c71b3e2SJacob Faibussowitsch   PetscCheck(M1 == M2 && N1 == N2, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_SIZ, "Non conforming matrix add: global sizes X %" PetscInt_FMT " x %" PetscInt_FMT ", Y %" PetscInt_FMT " x %" PetscInt_FMT, M1, N1, M2, N2);
1122c71b3e2SJacob Faibussowitsch   PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Non conforming matrix add: local sizes X %" PetscInt_FMT " x %" PetscInt_FMT ", Y %" PetscInt_FMT " x %" PetscInt_FMT, m1, n1, m2, n2);
1132c71b3e2SJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)");
1142c71b3e2SJacob Faibussowitsch   PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)");
1153ba16761SJacob Faibussowitsch   if (a == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
116edf9a46cSStefano Zampini   if (Y == X) {
1179566063dSJacob Faibussowitsch     PetscCall(MatScale(Y, 1.0 + a));
1183ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
119edf9a46cSStefano Zampini   }
120013e2dc7SBarry Smith   PetscCall(PetscObjectObjectTypeCompare((PetscObject)X, (PetscObject)Y, &sametype));
1219566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0));
12293c87e65SStefano Zampini   if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) {
123dbbe0bcdSBarry Smith     PetscUseTypeMethod(Y, axpy, a, X, str);
124d4bb536fSBarry Smith   } else {
125013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
126a683ea4aSPierre Jolivet     if (transpose) {
1279566063dSJacob Faibussowitsch       PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X));
128a683ea4aSPierre Jolivet     } else {
129013e2dc7SBarry Smith       PetscCall(PetscObjectTypeCompareAny((PetscObject)Y, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
130a683ea4aSPierre Jolivet       if (transpose) {
1319566063dSJacob Faibussowitsch         PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y));
132a683ea4aSPierre Jolivet       } else {
133*b20dac60SPierre Jolivet         PetscCall(MatAXPY_BasicWithTypeCompare(Y, a, X, str));
134a683ea4aSPierre Jolivet       }
135a683ea4aSPierre Jolivet     }
1362e37cc59SJose E. Roman   }
1379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0));
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139607cd303SBarry Smith }
140607cd303SBarry Smith 
141d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
142d71ae5a4SJacob Faibussowitsch {
14327afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL;
144646531bbSStefano Zampini 
145646531bbSStefano Zampini   PetscFunctionBegin;
146646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
1479566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall));
148646531bbSStefano Zampini   if (preall) {
1499566063dSJacob Faibussowitsch     PetscCall((*preall)(Y, X, B));
150646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
151646531bbSStefano Zampini     Mat      preallocator;
152646531bbSStefano Zampini     PetscInt r, rstart, rend;
153646531bbSStefano Zampini     PetscInt m, n, M, N;
154646531bbSStefano Zampini 
1559566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
1569566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
1579566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Y, &M, &N));
1589566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Y, &m, &n));
1599566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator));
1609566063dSJacob Faibussowitsch     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1619566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap));
1629566063dSJacob Faibussowitsch     PetscCall(MatSetUp(preallocator));
1639566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
164646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
165646531bbSStefano Zampini       PetscInt           ncols;
166646531bbSStefano Zampini       const PetscInt    *row;
167646531bbSStefano Zampini       const PetscScalar *vals;
168646531bbSStefano Zampini 
1699566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, r, &ncols, &row, &vals));
1709566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1719566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals));
1729566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, r, &ncols, &row, &vals));
1739566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1749566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals));
175646531bbSStefano Zampini     }
1769566063dSJacob Faibussowitsch     PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
1799566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
1809566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
181646531bbSStefano Zampini 
1829566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B));
1839566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name));
1849566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap));
1859566063dSJacob Faibussowitsch     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B));
1869566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&preallocator));
187646531bbSStefano Zampini   }
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
189646531bbSStefano Zampini }
190646531bbSStefano Zampini 
191d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str)
192d71ae5a4SJacob Faibussowitsch {
193edf9a46cSStefano Zampini   PetscFunctionBegin;
194*b20dac60SPierre Jolivet   if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) {
195*b20dac60SPierre Jolivet     PetscBool isdense;
196edf9a46cSStefano Zampini 
19793c87e65SStefano Zampini     /* no need to preallocate if Y is dense */
1989566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, ""));
199*b20dac60SPierre Jolivet     if (isdense) str = SUBSET_NONZERO_PATTERN;
200be705e3aSPierre Jolivet   }
201a6ab3590SBarry Smith   if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) {
202edf9a46cSStefano Zampini     PetscInt           i, start, end, j, ncols, m, n;
20338baddfdSBarry Smith     const PetscInt    *row;
204b3cc6726SBarry Smith     PetscScalar       *val;
205b3cc6726SBarry Smith     const PetscScalar *vals;
206607cd303SBarry Smith 
2079566063dSJacob Faibussowitsch     PetscCall(MatGetSize(X, &m, &n));
2089566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(X, &start, &end));
2099566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
210f4df32b1SMatthew Knepley     if (a == 1.0) {
211d4bb536fSBarry Smith       for (i = start; i < end; i++) {
2129566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2139566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES));
2149566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
215d4bb536fSBarry Smith       }
216d4bb536fSBarry Smith     } else {
21721a3365eSStefano Zampini       PetscInt vs = 100;
21821a3365eSStefano Zampini       /* realloc if needed, as this function may be used in parallel */
2199566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(vs, &val));
22006be10caSBarry Smith       for (i = start; i < end; i++) {
2219566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
22221a3365eSStefano Zampini         if (vs < ncols) {
22321a3365eSStefano Zampini           vs = PetscMin(2 * ncols, n);
2249566063dSJacob Faibussowitsch           PetscCall(PetscRealloc(vs * sizeof(*val), &val));
2256f79c3a4SBarry Smith         }
22621a3365eSStefano Zampini         for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2279566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES));
2289566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
2296f79c3a4SBarry Smith       }
2309566063dSJacob Faibussowitsch       PetscCall(PetscFree(val));
231d4bb536fSBarry Smith     }
2329566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
2339566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
2349566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
235edf9a46cSStefano Zampini   } else {
236edf9a46cSStefano Zampini     Mat B;
237edf9a46cSStefano Zampini 
2389566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B));
2399566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
2409566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y, &B));
241edf9a46cSStefano Zampini   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2436f79c3a4SBarry Smith }
244052efed2SBarry Smith 
245d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str)
246d71ae5a4SJacob Faibussowitsch {
247ec7775f6SShri Abhyankar   PetscInt           i, start, end, j, ncols, m, n;
248ec7775f6SShri Abhyankar   const PetscInt    *row;
249ec7775f6SShri Abhyankar   PetscScalar       *val;
250ec7775f6SShri Abhyankar   const PetscScalar *vals;
251ec7775f6SShri Abhyankar 
252ec7775f6SShri Abhyankar   PetscFunctionBegin;
2539566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &m, &n));
2549566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(X, &start, &end));
2559566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(Y));
2569566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(X));
257ec7775f6SShri Abhyankar   if (a == 1.0) {
258ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2599566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2609566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2619566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
262ec7775f6SShri Abhyankar 
2639566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2649566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2659566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
266ec7775f6SShri Abhyankar     }
267ec7775f6SShri Abhyankar   } else {
26821a3365eSStefano Zampini     PetscInt vs = 100;
26921a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
2709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(vs, &val));
271ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2729566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2739566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2749566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
275ec7775f6SShri Abhyankar 
2769566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
27721a3365eSStefano Zampini       if (vs < ncols) {
27821a3365eSStefano Zampini         vs = PetscMin(2 * ncols, n);
2799566063dSJacob Faibussowitsch         PetscCall(PetscRealloc(vs * sizeof(*val), &val));
280ec7775f6SShri Abhyankar       }
28121a3365eSStefano Zampini       for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2829566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES));
2839566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
284ec7775f6SShri Abhyankar     }
2859566063dSJacob Faibussowitsch     PetscCall(PetscFree(val));
286ec7775f6SShri Abhyankar   }
2879566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(Y));
2889566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(X));
2899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
292ec7775f6SShri Abhyankar }
293ec7775f6SShri Abhyankar 
294052efed2SBarry Smith /*@
2952ef1f0ffSBarry Smith   MatShift - Computes `Y =  Y + a I`, where `a` is a `PetscScalar`
296052efed2SBarry Smith 
297c3339decSBarry Smith   Neighbor-wise Collective
298fee21e36SBarry Smith 
29998a79cdbSBarry Smith   Input Parameters:
30027430b45SBarry Smith + Y - the matrix
30127430b45SBarry Smith - a - the `PetscScalar`
30298a79cdbSBarry Smith 
3032860a424SLois Curfman McInnes   Level: intermediate
3042860a424SLois Curfman McInnes 
30595452b02SPatrick Sanan   Notes:
30627430b45SBarry Smith   If `Y` is a rectangular matrix, the shift is done on the main diagonal of the matrix (https://en.wikipedia.org/wiki/Main_diagonal)
30784e19e3eSJunchao Zhang 
3082ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3096f33a894SBarry Smith   fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an
310edf9a46cSStefano Zampini   entry). No operation is performed when a is zero.
3116f33a894SBarry Smith 
31227430b45SBarry Smith   To form Y = Y + diag(V) use `MatDiagonalSet()`
3130c0fd78eSBarry Smith 
3141cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()`
315052efed2SBarry Smith  @*/
316d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift(Mat Y, PetscScalar a)
317d71ae5a4SJacob Faibussowitsch {
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3190700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3205f80ce2aSJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3215f80ce2aSJacob Faibussowitsch   PetscCheck(!Y->factortype, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3224994cf47SJed Brown   MatCheckPreallocated(Y, 1);
3233ba16761SJacob Faibussowitsch   if (a == 0.0) PetscFunctionReturn(PETSC_SUCCESS);
324b50b34bdSBarry Smith 
325dbbe0bcdSBarry Smith   if (Y->ops->shift) PetscUseTypeMethod(Y, shift, a);
3269566063dSJacob Faibussowitsch   else PetscCall(MatShift_Basic(Y, a));
3277d68702bSBarry Smith 
3289566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
330052efed2SBarry Smith }
3316d84be18SBarry Smith 
332d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is)
333d71ae5a4SJacob Faibussowitsch {
33467576b8bSBarry Smith   PetscInt           i, start, end;
335fa2e578bSStefano Zampini   const PetscScalar *v;
33609f38230SBarry Smith 
33709f38230SBarry Smith   PetscFunctionBegin;
3389566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y, &start, &end));
3399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(D, &v));
34048a46eb9SPierre Jolivet   for (i = start; i < end; i++) PetscCall(MatSetValues(Y, 1, &i, 1, &i, v + i - start, is));
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(D, &v));
3429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
3439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34509f38230SBarry Smith }
34609f38230SBarry Smith 
3476d84be18SBarry Smith /*@
3482ef1f0ffSBarry Smith   MatDiagonalSet - Computes `Y` = `Y` + `D`, where `D` is a diagonal matrix
3492ef1f0ffSBarry Smith   that is represented as a vector. Or Y[i,i] = D[i] if `InsertMode` is
3502ef1f0ffSBarry Smith   `INSERT_VALUES`.
3516d84be18SBarry Smith 
352c3339decSBarry Smith   Neighbor-wise Collective
35348e586daSJose E. Roman 
3546d84be18SBarry Smith   Input Parameters:
35598a79cdbSBarry Smith + Y  - the input matrix
356f56f2b3fSBarry Smith . D  - the diagonal matrix, represented as a vector
357fe59aa6dSJacob Faibussowitsch - is - `INSERT_VALUES` or `ADD_VALUES`
3586f33a894SBarry Smith 
3592860a424SLois Curfman McInnes   Level: intermediate
3602860a424SLois Curfman McInnes 
3612ef1f0ffSBarry Smith   Note:
3622ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3632ef1f0ffSBarry Smith   fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an
3642ef1f0ffSBarry Smith   entry).
3652ef1f0ffSBarry Smith 
3661cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatShift()`, `MatScale()`, `MatDiagonalScale()`
3676d84be18SBarry Smith @*/
368d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is)
369d71ae5a4SJacob Faibussowitsch {
37067576b8bSBarry Smith   PetscInt matlocal, veclocal;
3716d84be18SBarry Smith 
3723a40ed3dSBarry Smith   PetscFunctionBegin;
3730700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3740700a824SBarry Smith   PetscValidHeaderSpecific(D, VEC_CLASSID, 2);
3752fbb122aSJose E. Roman   MatCheckPreallocated(Y, 1);
3769566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &matlocal, NULL));
3779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(D, &veclocal));
3785f80ce2aSJacob Faibussowitsch   PetscCheck(matlocal == veclocal, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number local rows of matrix %" PetscInt_FMT " does not match that of vector for diagonal %" PetscInt_FMT, matlocal, veclocal);
379dbbe0bcdSBarry Smith   if (Y->ops->diagonalset) PetscUseTypeMethod(Y, diagonalset, D, is);
380dbbe0bcdSBarry Smith   else PetscCall(MatDiagonalSet_Default(Y, D, is));
3819566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3836d84be18SBarry Smith }
384d4bb536fSBarry Smith 
385d4bb536fSBarry Smith /*@
38604aac2b0SHong Zhang   MatAYPX - Computes Y = a*Y + X.
387d4bb536fSBarry Smith 
3882ef1f0ffSBarry Smith   Logically Collective
389fee21e36SBarry Smith 
39098a79cdbSBarry Smith   Input Parameters:
3912ef1f0ffSBarry Smith + a   - the `PetscScalar` multiplier
39204aac2b0SHong Zhang . Y   - the first matrix
39304aac2b0SHong Zhang . X   - the second matrix
394aaa8cc7dSPierre Jolivet - str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN`, `UNKNOWN_NONZERO_PATTERN`, or `SUBSET_NONZERO_PATTERN` (nonzeros of `X` is a subset of `Y`'s)
395d4bb536fSBarry Smith 
3962860a424SLois Curfman McInnes   Level: intermediate
3972860a424SLois Curfman McInnes 
3981cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAXPY()`
399d4bb536fSBarry Smith  @*/
400d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAYPX(Mat Y, PetscScalar a, Mat X, MatStructure str)
401d71ae5a4SJacob Faibussowitsch {
4023a40ed3dSBarry Smith   PetscFunctionBegin;
4039566063dSJacob Faibussowitsch   PetscCall(MatScale(Y, a));
4049566063dSJacob Faibussowitsch   PetscCall(MatAXPY(Y, 1.0, X, str));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406d4bb536fSBarry Smith }
407b0a32e0cSBarry Smith 
4085d83a8b1SBarry Smith /*@
4090bacdadaSStefano Zampini   MatComputeOperator - Computes the explicit matrix
410b0a32e0cSBarry Smith 
411c3339decSBarry Smith   Collective
412b0a32e0cSBarry Smith 
413d8d19677SJose E. Roman   Input Parameters:
414186a3e20SStefano Zampini + inmat   - the matrix
415186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator
416b0a32e0cSBarry Smith 
417b0a32e0cSBarry Smith   Output Parameter:
418a5b23f4aSJose E. Roman . mat - the explicit  operator
419b0a32e0cSBarry Smith 
4202ef1f0ffSBarry Smith   Level: advanced
4212ef1f0ffSBarry Smith 
42211a5261eSBarry Smith   Note:
423f13dfd9eSBarry Smith   This computation is done by applying the operator to columns of the identity matrix.
424186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
4252ef1f0ffSBarry Smith   Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
426b0a32e0cSBarry Smith 
4271cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperatorTranspose()`
428b0a32e0cSBarry Smith @*/
429d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperator(Mat inmat, MatType mattype, Mat *mat)
430d71ae5a4SJacob Faibussowitsch {
431b0a32e0cSBarry Smith   PetscFunctionBegin;
4320700a824SBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
4334f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
4349566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(inmat, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43624f910e3SHong Zhang }
4374325cce7SMatthew G Knepley 
4385d83a8b1SBarry Smith /*@
4390bacdadaSStefano Zampini   MatComputeOperatorTranspose - Computes the explicit matrix representation of
4402ef1f0ffSBarry Smith   a give matrix that can apply `MatMultTranspose()`
441f3b1f45cSBarry Smith 
442c3339decSBarry Smith   Collective
443f3b1f45cSBarry Smith 
4446b867d5aSJose E. Roman   Input Parameters:
4456b867d5aSJose E. Roman + inmat   - the matrix
4466b867d5aSJose E. Roman - mattype - the matrix type for the explicit operator
447f3b1f45cSBarry Smith 
448f3b1f45cSBarry Smith   Output Parameter:
449a5b23f4aSJose E. Roman . mat - the explicit  operator transposed
450f3b1f45cSBarry Smith 
4512ef1f0ffSBarry Smith   Level: advanced
4522ef1f0ffSBarry Smith 
45311a5261eSBarry Smith   Note:
454186a3e20SStefano Zampini   This computation is done by applying the transpose of the operator to columns of the identity matrix.
455186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
4562ef1f0ffSBarry Smith   Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
457f3b1f45cSBarry Smith 
4581cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperator()`
459f3b1f45cSBarry Smith @*/
460d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperatorTranspose(Mat inmat, MatType mattype, Mat *mat)
461d71ae5a4SJacob Faibussowitsch {
462186a3e20SStefano Zampini   Mat A;
463f3b1f45cSBarry Smith 
464f3b1f45cSBarry Smith   PetscFunctionBegin;
465f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
4664f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
4679566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(inmat, &A));
4689566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(A, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471f3b1f45cSBarry Smith }
472f3b1f45cSBarry Smith 
473f3b1f45cSBarry Smith /*@
4742ce66baaSPierre Jolivet   MatFilter - Set all values in the matrix with an absolute value less than or equal to the tolerance to zero, and optionally compress the underlying storage
4754325cce7SMatthew G Knepley 
4764325cce7SMatthew G Knepley   Input Parameters:
4774325cce7SMatthew G Knepley + A        - The matrix
4782ce66baaSPierre Jolivet . tol      - The zero tolerance
4792ce66baaSPierre Jolivet . compress - Whether the storage from the input matrix `A` should be compressed once values less than or equal to `tol` are set to zero
4802ce66baaSPierre Jolivet - keep     - If `compress` is true and for a given row of `A`, the diagonal coefficient is less than or equal to `tol`, indicates whether it should be left in the structure or eliminated as well
4814325cce7SMatthew G Knepley 
4824325cce7SMatthew G Knepley   Level: intermediate
4834325cce7SMatthew G Knepley 
484d9b70fcdSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatZeroEntries()`, `MatEliminateZeros()`, `VecFilter()`
4853fc99919SSatish Balay  @*/
4862ce66baaSPierre Jolivet PetscErrorCode MatFilter(Mat A, PetscReal tol, PetscBool compress, PetscBool keep)
487d71ae5a4SJacob Faibussowitsch {
488038df967SPierre Jolivet   Mat          a;
4894325cce7SMatthew G Knepley   PetscScalar *newVals;
49034bf4ff8Smarkadams4   PetscInt    *newCols, rStart, rEnd, maxRows, r, colMax = 0, nnz0 = 0, nnz1 = 0;
491038df967SPierre Jolivet   PetscBool    flg;
4924325cce7SMatthew G Knepley 
4934325cce7SMatthew G Knepley   PetscFunctionBegin;
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, ""));
495038df967SPierre Jolivet   if (flg) {
4969566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLocalMatrix(A, &a));
4979566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(a, &r));
4989566063dSJacob Faibussowitsch     PetscCall(MatGetSize(a, &rStart, &rEnd));
4999566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(a, &newVals));
500038df967SPierre Jolivet     for (; colMax < rEnd; ++colMax) {
50182fa6941SPierre Jolivet       for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) <= tol ? 0.0 : newVals[maxRows + colMax * r];
502038df967SPierre Jolivet     }
5039566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(a, &newVals));
504038df967SPierre Jolivet   } else {
5053da7c217SPierre Jolivet     const PetscInt *ranges;
5063da7c217SPierre Jolivet     PetscMPIInt     rank, size;
5073da7c217SPierre Jolivet 
5083da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
5093da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
5103da7c217SPierre Jolivet     PetscCall(MatGetOwnershipRanges(A, &ranges));
5113da7c217SPierre Jolivet     rStart = ranges[rank];
5123da7c217SPierre Jolivet     rEnd   = ranges[rank + 1];
5139566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
5144325cce7SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
5154325cce7SMatthew G Knepley       PetscInt ncols;
5164325cce7SMatthew G Knepley 
5179566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, r, &ncols, NULL, NULL));
5185f80ce2aSJacob Faibussowitsch       colMax = PetscMax(colMax, ncols);
5199566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL));
5204325cce7SMatthew G Knepley     }
5213da7c217SPierre Jolivet     maxRows = 0;
5223da7c217SPierre Jolivet     for (r = 0; r < size; ++r) maxRows = PetscMax(maxRows, ranges[r + 1] - ranges[r]);
5233da7c217SPierre Jolivet     PetscCall(PetscCalloc2(colMax, &newCols, colMax, &newVals));
5249566063dSJacob Faibussowitsch     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */
5259566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
526cf5d3338SPierre Jolivet     /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd()             */
527cf5d3338SPierre Jolivet     /* that are potentially called many times depending on the distribution of A */
5284325cce7SMatthew G Knepley     for (r = rStart; r < rStart + maxRows; ++r) {
5293da7c217SPierre Jolivet       if (r < rEnd) {
5304325cce7SMatthew G Knepley         const PetscScalar *vals;
5314325cce7SMatthew G Knepley         const PetscInt    *cols;
5323da7c217SPierre Jolivet         PetscInt           ncols, newcols = 0, c;
5334325cce7SMatthew G Knepley 
5349566063dSJacob Faibussowitsch         PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
53534bf4ff8Smarkadams4         nnz0 += ncols - 1;
5364325cce7SMatthew G Knepley         for (c = 0; c < ncols; ++c) {
53782fa6941SPierre Jolivet           if (PetscUnlikely(PetscAbsScalar(vals[c]) <= tol)) newCols[newcols++] = cols[c];
5384325cce7SMatthew G Knepley         }
53934bf4ff8Smarkadams4         nnz1 += ncols - newcols - 1;
5409566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
5419566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
5424325cce7SMatthew G Knepley       }
5439566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5449566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5454325cce7SMatthew G Knepley     }
5469566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
5479566063dSJacob Faibussowitsch     PetscCall(PetscFree2(newCols, newVals));
5489566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */
54934bf4ff8Smarkadams4     if (nnz0 > 0) PetscCall(PetscInfo(NULL, "Filtering left %g%% edges in graph\n", 100 * (double)nnz1 / (double)nnz0));
5505afa97f7SPierre Jolivet     else PetscCall(PetscInfo(NULL, "Warning: %" PetscInt_FMT " edges to filter with %" PetscInt_FMT " rows\n", nnz0, maxRows));
551038df967SPierre Jolivet   }
5522ce66baaSPierre Jolivet   if (compress && A->ops->eliminatezeros) {
5532ce66baaSPierre Jolivet     Mat       B;
55424c7e8a4SPierre Jolivet     PetscBool flg;
5552ce66baaSPierre Jolivet 
55624c7e8a4SPierre Jolivet     PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJHIPSPARSE, MATMPIAIJHIPSPARSE, ""));
55724c7e8a4SPierre Jolivet     if (!flg) {
5582ce66baaSPierre Jolivet       PetscCall(MatEliminateZeros(A, keep));
5592ce66baaSPierre Jolivet       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
5602ce66baaSPierre Jolivet       PetscCall(MatHeaderReplace(A, &B));
5612ce66baaSPierre Jolivet     }
56224c7e8a4SPierre Jolivet   }
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5644325cce7SMatthew G Knepley }
565