xref: /petsc/src/mat/utils/axpy.c (revision 4166f22e6d0bce267963a27c503b606fd4ea00f5)
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;
6*4166f22eSPierre Jolivet   PetscScalar vshift, vscale;
75f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat, Mat *);
8d54c938cSPierre Jolivet 
9d54c938cSPierre Jolivet   PetscFunctionBegin;
10*4166f22eSPierre 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));
11*4166f22eSPierre Jolivet   else {
12*4166f22eSPierre Jolivet     vshift = 0.0;
13*4166f22eSPierre Jolivet     vscale = 1.0;
14*4166f22eSPierre 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   }
37*4166f22eSPierre Jolivet   PetscCall(MatAXPY(A, a * vscale, F, str));
38*4166f22eSPierre Jolivet   PetscCall(MatShift(A, a * vshift));
399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41d54c938cSPierre Jolivet }
42d54c938cSPierre Jolivet 
4306be10caSBarry Smith /*@
4421c89e3eSBarry Smith   MatAXPY - Computes Y = a*X + Y.
456f79c3a4SBarry Smith 
46c3339decSBarry Smith   Logically Collective
47fee21e36SBarry Smith 
4898a79cdbSBarry Smith   Input Parameters:
49607cd303SBarry Smith + a   - the scalar multiplier
50607cd303SBarry Smith . X   - the first matrix
51607cd303SBarry Smith . Y   - the second matrix
5227430b45SBarry 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)
53edf9a46cSStefano Zampini 
542860a424SLois Curfman McInnes   Level: intermediate
552860a424SLois Curfman McInnes 
561cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAYPX()`
5706be10caSBarry Smith  @*/
58d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str)
59d71ae5a4SJacob Faibussowitsch {
60646531bbSStefano Zampini   PetscInt  M1, M2, N1, N2;
61c1ac3661SBarry Smith   PetscInt  m1, m2, n1, n2;
62a683ea4aSPierre Jolivet   PetscBool sametype, transpose;
636f79c3a4SBarry Smith 
643a40ed3dSBarry Smith   PetscFunctionBegin;
650700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
66c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y, a, 2);
67646531bbSStefano Zampini   PetscValidHeaderSpecific(X, MAT_CLASSID, 3);
68646531bbSStefano Zampini   PetscCheckSameComm(Y, 1, X, 3);
699566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M1, &N1));
709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M2, &N2));
719566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m1, &n1));
729566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m2, &n2));
732c71b3e2SJacob 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);
742c71b3e2SJacob 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);
752c71b3e2SJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)");
762c71b3e2SJacob Faibussowitsch   PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)");
773ba16761SJacob Faibussowitsch   if (a == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
78edf9a46cSStefano Zampini   if (Y == X) {
799566063dSJacob Faibussowitsch     PetscCall(MatScale(Y, 1.0 + a));
803ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
81edf9a46cSStefano Zampini   }
82013e2dc7SBarry Smith   PetscCall(PetscObjectObjectTypeCompare((PetscObject)X, (PetscObject)Y, &sametype));
839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0));
8493c87e65SStefano Zampini   if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) {
85dbbe0bcdSBarry Smith     PetscUseTypeMethod(Y, axpy, a, X, str);
86d4bb536fSBarry Smith   } else {
87013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
88a683ea4aSPierre Jolivet     if (transpose) {
899566063dSJacob Faibussowitsch       PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X));
90a683ea4aSPierre Jolivet     } else {
91013e2dc7SBarry Smith       PetscCall(PetscObjectTypeCompareAny((PetscObject)Y, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
92a683ea4aSPierre Jolivet       if (transpose) {
939566063dSJacob Faibussowitsch         PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y));
94a683ea4aSPierre Jolivet       } else {
959566063dSJacob Faibussowitsch         PetscCall(MatAXPY_Basic(Y, a, X, str));
96607cd303SBarry Smith       }
97a683ea4aSPierre Jolivet     }
98a683ea4aSPierre Jolivet   }
999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
101607cd303SBarry Smith }
102607cd303SBarry Smith 
103d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
104d71ae5a4SJacob Faibussowitsch {
10527afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL;
106646531bbSStefano Zampini 
107646531bbSStefano Zampini   PetscFunctionBegin;
108646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall));
110646531bbSStefano Zampini   if (preall) {
1119566063dSJacob Faibussowitsch     PetscCall((*preall)(Y, X, B));
112646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
113646531bbSStefano Zampini     Mat      preallocator;
114646531bbSStefano Zampini     PetscInt r, rstart, rend;
115646531bbSStefano Zampini     PetscInt m, n, M, N;
116646531bbSStefano Zampini 
1179566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
1189566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
1199566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Y, &M, &N));
1209566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Y, &m, &n));
1219566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator));
1229566063dSJacob Faibussowitsch     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1239566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap));
1249566063dSJacob Faibussowitsch     PetscCall(MatSetUp(preallocator));
1259566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
126646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
127646531bbSStefano Zampini       PetscInt           ncols;
128646531bbSStefano Zampini       const PetscInt    *row;
129646531bbSStefano Zampini       const PetscScalar *vals;
130646531bbSStefano Zampini 
1319566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, r, &ncols, &row, &vals));
1329566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1339566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals));
1349566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, r, &ncols, &row, &vals));
1359566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1369566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals));
137646531bbSStefano Zampini     }
1389566063dSJacob Faibussowitsch     PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1399566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1409566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
1419566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
1429566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
143646531bbSStefano Zampini 
1449566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B));
1459566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name));
1469566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap));
1479566063dSJacob Faibussowitsch     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B));
1489566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&preallocator));
149646531bbSStefano Zampini   }
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151646531bbSStefano Zampini }
152646531bbSStefano Zampini 
153d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str)
154d71ae5a4SJacob Faibussowitsch {
155be705e3aSPierre Jolivet   PetscBool isshell, isdense, isnest;
156edf9a46cSStefano Zampini 
157edf9a46cSStefano Zampini   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   PetscCall(MatIsShell(Y, &isshell));
159edf9a46cSStefano Zampini   if (isshell) { /* MatShell has special support for AXPY */
160edf9a46cSStefano Zampini     PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure);
161edf9a46cSStefano Zampini 
1629566063dSJacob Faibussowitsch     PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f));
163edf9a46cSStefano Zampini     if (f) {
1649566063dSJacob Faibussowitsch       PetscCall((*f)(Y, a, X, str));
1653ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
166edf9a46cSStefano Zampini     }
167edf9a46cSStefano Zampini   }
16893c87e65SStefano Zampini   /* no need to preallocate if Y is dense */
1699566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, ""));
170be705e3aSPierre Jolivet   if (isdense) {
1719566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &isnest));
172be705e3aSPierre Jolivet     if (isnest) {
1739566063dSJacob Faibussowitsch       PetscCall(MatAXPY_Dense_Nest(Y, a, X));
1743ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
175be705e3aSPierre Jolivet     }
176be705e3aSPierre Jolivet     if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN;
177be705e3aSPierre Jolivet   }
178a6ab3590SBarry Smith   if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) {
179edf9a46cSStefano Zampini     PetscInt           i, start, end, j, ncols, m, n;
18038baddfdSBarry Smith     const PetscInt    *row;
181b3cc6726SBarry Smith     PetscScalar       *val;
182b3cc6726SBarry Smith     const PetscScalar *vals;
183607cd303SBarry Smith 
1849566063dSJacob Faibussowitsch     PetscCall(MatGetSize(X, &m, &n));
1859566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(X, &start, &end));
1869566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
187f4df32b1SMatthew Knepley     if (a == 1.0) {
188d4bb536fSBarry Smith       for (i = start; i < end; i++) {
1899566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
1909566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES));
1919566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
192d4bb536fSBarry Smith       }
193d4bb536fSBarry Smith     } else {
19421a3365eSStefano Zampini       PetscInt vs = 100;
19521a3365eSStefano Zampini       /* realloc if needed, as this function may be used in parallel */
1969566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(vs, &val));
19706be10caSBarry Smith       for (i = start; i < end; i++) {
1989566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
19921a3365eSStefano Zampini         if (vs < ncols) {
20021a3365eSStefano Zampini           vs = PetscMin(2 * ncols, n);
2019566063dSJacob Faibussowitsch           PetscCall(PetscRealloc(vs * sizeof(*val), &val));
2026f79c3a4SBarry Smith         }
20321a3365eSStefano Zampini         for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2049566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES));
2059566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
2066f79c3a4SBarry Smith       }
2079566063dSJacob Faibussowitsch       PetscCall(PetscFree(val));
208d4bb536fSBarry Smith     }
2099566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
2109566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
2119566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
212edf9a46cSStefano Zampini   } else {
213edf9a46cSStefano Zampini     Mat B;
214edf9a46cSStefano Zampini 
2159566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B));
2169566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
2179566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y, &B));
218edf9a46cSStefano Zampini   }
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2206f79c3a4SBarry Smith }
221052efed2SBarry Smith 
222d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str)
223d71ae5a4SJacob Faibussowitsch {
224ec7775f6SShri Abhyankar   PetscInt           i, start, end, j, ncols, m, n;
225ec7775f6SShri Abhyankar   const PetscInt    *row;
226ec7775f6SShri Abhyankar   PetscScalar       *val;
227ec7775f6SShri Abhyankar   const PetscScalar *vals;
228ec7775f6SShri Abhyankar 
229ec7775f6SShri Abhyankar   PetscFunctionBegin;
2309566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &m, &n));
2319566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(X, &start, &end));
2329566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(Y));
2339566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(X));
234ec7775f6SShri Abhyankar   if (a == 1.0) {
235ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2369566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2379566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2389566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
239ec7775f6SShri Abhyankar 
2409566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2419566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2429566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
243ec7775f6SShri Abhyankar     }
244ec7775f6SShri Abhyankar   } else {
24521a3365eSStefano Zampini     PetscInt vs = 100;
24621a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
2479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(vs, &val));
248ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2499566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2509566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2519566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
252ec7775f6SShri Abhyankar 
2539566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
25421a3365eSStefano Zampini       if (vs < ncols) {
25521a3365eSStefano Zampini         vs = PetscMin(2 * ncols, n);
2569566063dSJacob Faibussowitsch         PetscCall(PetscRealloc(vs * sizeof(*val), &val));
257ec7775f6SShri Abhyankar       }
25821a3365eSStefano Zampini       for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2599566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES));
2609566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
261ec7775f6SShri Abhyankar     }
2629566063dSJacob Faibussowitsch     PetscCall(PetscFree(val));
263ec7775f6SShri Abhyankar   }
2649566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(Y));
2659566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(X));
2669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
269ec7775f6SShri Abhyankar }
270ec7775f6SShri Abhyankar 
271052efed2SBarry Smith /*@
2722ef1f0ffSBarry Smith   MatShift - Computes `Y =  Y + a I`, where `a` is a `PetscScalar`
273052efed2SBarry Smith 
274c3339decSBarry Smith   Neighbor-wise Collective
275fee21e36SBarry Smith 
27698a79cdbSBarry Smith   Input Parameters:
27727430b45SBarry Smith + Y - the matrix
27827430b45SBarry Smith - a - the `PetscScalar`
27998a79cdbSBarry Smith 
2802860a424SLois Curfman McInnes   Level: intermediate
2812860a424SLois Curfman McInnes 
28295452b02SPatrick Sanan   Notes:
28327430b45SBarry 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)
28484e19e3eSJunchao Zhang 
2852ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2866f33a894SBarry 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
287edf9a46cSStefano Zampini   entry). No operation is performed when a is zero.
2886f33a894SBarry Smith 
28927430b45SBarry Smith   To form Y = Y + diag(V) use `MatDiagonalSet()`
2900c0fd78eSBarry Smith 
2911cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()`
292052efed2SBarry Smith  @*/
293d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift(Mat Y, PetscScalar a)
294d71ae5a4SJacob Faibussowitsch {
2953a40ed3dSBarry Smith   PetscFunctionBegin;
2960700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
2975f80ce2aSJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2985f80ce2aSJacob Faibussowitsch   PetscCheck(!Y->factortype, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2994994cf47SJed Brown   MatCheckPreallocated(Y, 1);
3003ba16761SJacob Faibussowitsch   if (a == 0.0) PetscFunctionReturn(PETSC_SUCCESS);
301b50b34bdSBarry Smith 
302dbbe0bcdSBarry Smith   if (Y->ops->shift) PetscUseTypeMethod(Y, shift, a);
3039566063dSJacob Faibussowitsch   else PetscCall(MatShift_Basic(Y, a));
3047d68702bSBarry Smith 
3059566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307052efed2SBarry Smith }
3086d84be18SBarry Smith 
309d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is)
310d71ae5a4SJacob Faibussowitsch {
31167576b8bSBarry Smith   PetscInt           i, start, end;
312fa2e578bSStefano Zampini   const PetscScalar *v;
31309f38230SBarry Smith 
31409f38230SBarry Smith   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y, &start, &end));
3169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(D, &v));
31748a46eb9SPierre Jolivet   for (i = start; i < end; i++) PetscCall(MatSetValues(Y, 1, &i, 1, &i, v + i - start, is));
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(D, &v));
3199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
3209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32209f38230SBarry Smith }
32309f38230SBarry Smith 
3246d84be18SBarry Smith /*@
3252ef1f0ffSBarry Smith   MatDiagonalSet - Computes `Y` = `Y` + `D`, where `D` is a diagonal matrix
3262ef1f0ffSBarry Smith   that is represented as a vector. Or Y[i,i] = D[i] if `InsertMode` is
3272ef1f0ffSBarry Smith   `INSERT_VALUES`.
3286d84be18SBarry Smith 
329c3339decSBarry Smith   Neighbor-wise Collective
33048e586daSJose E. Roman 
3316d84be18SBarry Smith   Input Parameters:
33298a79cdbSBarry Smith + Y  - the input matrix
333f56f2b3fSBarry Smith . D  - the diagonal matrix, represented as a vector
334fe59aa6dSJacob Faibussowitsch - is - `INSERT_VALUES` or `ADD_VALUES`
3356f33a894SBarry Smith 
3362860a424SLois Curfman McInnes   Level: intermediate
3372860a424SLois Curfman McInnes 
3382ef1f0ffSBarry Smith   Note:
3392ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3402ef1f0ffSBarry 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
3412ef1f0ffSBarry Smith   entry).
3422ef1f0ffSBarry Smith 
3431cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatShift()`, `MatScale()`, `MatDiagonalScale()`
3446d84be18SBarry Smith @*/
345d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is)
346d71ae5a4SJacob Faibussowitsch {
34767576b8bSBarry Smith   PetscInt matlocal, veclocal;
3486d84be18SBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
3500700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3510700a824SBarry Smith   PetscValidHeaderSpecific(D, VEC_CLASSID, 2);
3522fbb122aSJose E. Roman   MatCheckPreallocated(Y, 1);
3539566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &matlocal, NULL));
3549566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(D, &veclocal));
3555f80ce2aSJacob 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);
356dbbe0bcdSBarry Smith   if (Y->ops->diagonalset) PetscUseTypeMethod(Y, diagonalset, D, is);
357dbbe0bcdSBarry Smith   else PetscCall(MatDiagonalSet_Default(Y, D, is));
3589566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3606d84be18SBarry Smith }
361d4bb536fSBarry Smith 
362d4bb536fSBarry Smith /*@
36304aac2b0SHong Zhang   MatAYPX - Computes Y = a*Y + X.
364d4bb536fSBarry Smith 
3652ef1f0ffSBarry Smith   Logically Collective
366fee21e36SBarry Smith 
36798a79cdbSBarry Smith   Input Parameters:
3682ef1f0ffSBarry Smith + a   - the `PetscScalar` multiplier
36904aac2b0SHong Zhang . Y   - the first matrix
37004aac2b0SHong Zhang . X   - the second matrix
371aaa8cc7dSPierre 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)
372d4bb536fSBarry Smith 
3732860a424SLois Curfman McInnes   Level: intermediate
3742860a424SLois Curfman McInnes 
3751cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAXPY()`
376d4bb536fSBarry Smith  @*/
377d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAYPX(Mat Y, PetscScalar a, Mat X, MatStructure str)
378d71ae5a4SJacob Faibussowitsch {
3793a40ed3dSBarry Smith   PetscFunctionBegin;
3809566063dSJacob Faibussowitsch   PetscCall(MatScale(Y, a));
3819566063dSJacob Faibussowitsch   PetscCall(MatAXPY(Y, 1.0, X, str));
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383d4bb536fSBarry Smith }
384b0a32e0cSBarry Smith 
385fe59aa6dSJacob Faibussowitsch /*@C
3860bacdadaSStefano Zampini   MatComputeOperator - Computes the explicit matrix
387b0a32e0cSBarry Smith 
388c3339decSBarry Smith   Collective
389b0a32e0cSBarry Smith 
390d8d19677SJose E. Roman   Input Parameters:
391186a3e20SStefano Zampini + inmat   - the matrix
392186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator
393b0a32e0cSBarry Smith 
394b0a32e0cSBarry Smith   Output Parameter:
395a5b23f4aSJose E. Roman . mat - the explicit  operator
396b0a32e0cSBarry Smith 
3972ef1f0ffSBarry Smith   Level: advanced
3982ef1f0ffSBarry Smith 
39911a5261eSBarry Smith   Note:
400f13dfd9eSBarry Smith   This computation is done by applying the operator to columns of the identity matrix.
401186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
4022ef1f0ffSBarry Smith   Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
403b0a32e0cSBarry Smith 
4041cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperatorTranspose()`
405b0a32e0cSBarry Smith @*/
406d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperator(Mat inmat, MatType mattype, Mat *mat)
407d71ae5a4SJacob Faibussowitsch {
408b0a32e0cSBarry Smith   PetscFunctionBegin;
4090700a824SBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
4104f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
4119566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(inmat, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41324f910e3SHong Zhang }
4144325cce7SMatthew G Knepley 
415fe59aa6dSJacob Faibussowitsch /*@C
4160bacdadaSStefano Zampini   MatComputeOperatorTranspose - Computes the explicit matrix representation of
4172ef1f0ffSBarry Smith   a give matrix that can apply `MatMultTranspose()`
418f3b1f45cSBarry Smith 
419c3339decSBarry Smith   Collective
420f3b1f45cSBarry Smith 
4216b867d5aSJose E. Roman   Input Parameters:
4226b867d5aSJose E. Roman + inmat   - the matrix
4236b867d5aSJose E. Roman - mattype - the matrix type for the explicit operator
424f3b1f45cSBarry Smith 
425f3b1f45cSBarry Smith   Output Parameter:
426a5b23f4aSJose E. Roman . mat - the explicit  operator transposed
427f3b1f45cSBarry Smith 
4282ef1f0ffSBarry Smith   Level: advanced
4292ef1f0ffSBarry Smith 
43011a5261eSBarry Smith   Note:
431186a3e20SStefano Zampini   This computation is done by applying the transpose of the operator to columns of the identity matrix.
432186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
4332ef1f0ffSBarry Smith   Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
434f3b1f45cSBarry Smith 
4351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperator()`
436f3b1f45cSBarry Smith @*/
437d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperatorTranspose(Mat inmat, MatType mattype, Mat *mat)
438d71ae5a4SJacob Faibussowitsch {
439186a3e20SStefano Zampini   Mat A;
440f3b1f45cSBarry Smith 
441f3b1f45cSBarry Smith   PetscFunctionBegin;
442f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
4434f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
4449566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(inmat, &A));
4459566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(A, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
448f3b1f45cSBarry Smith }
449f3b1f45cSBarry Smith 
450f3b1f45cSBarry Smith /*@
4512ce66baaSPierre 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
4524325cce7SMatthew G Knepley 
4534325cce7SMatthew G Knepley   Input Parameters:
4544325cce7SMatthew G Knepley + A        - The matrix
4552ce66baaSPierre Jolivet . tol      - The zero tolerance
4562ce66baaSPierre 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
4572ce66baaSPierre 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
4584325cce7SMatthew G Knepley 
4594325cce7SMatthew G Knepley   Level: intermediate
4604325cce7SMatthew G Knepley 
461d9b70fcdSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatZeroEntries()`, `MatEliminateZeros()`, `VecFilter()`
4623fc99919SSatish Balay  @*/
4632ce66baaSPierre Jolivet PetscErrorCode MatFilter(Mat A, PetscReal tol, PetscBool compress, PetscBool keep)
464d71ae5a4SJacob Faibussowitsch {
465038df967SPierre Jolivet   Mat          a;
4664325cce7SMatthew G Knepley   PetscScalar *newVals;
46734bf4ff8Smarkadams4   PetscInt    *newCols, rStart, rEnd, maxRows, r, colMax = 0, nnz0 = 0, nnz1 = 0;
468038df967SPierre Jolivet   PetscBool    flg;
4694325cce7SMatthew G Knepley 
4704325cce7SMatthew G Knepley   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, ""));
472038df967SPierre Jolivet   if (flg) {
4739566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLocalMatrix(A, &a));
4749566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(a, &r));
4759566063dSJacob Faibussowitsch     PetscCall(MatGetSize(a, &rStart, &rEnd));
4769566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(a, &newVals));
477038df967SPierre Jolivet     for (; colMax < rEnd; ++colMax) {
47882fa6941SPierre Jolivet       for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) <= tol ? 0.0 : newVals[maxRows + colMax * r];
479038df967SPierre Jolivet     }
4809566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(a, &newVals));
481038df967SPierre Jolivet   } else {
4823da7c217SPierre Jolivet     const PetscInt *ranges;
4833da7c217SPierre Jolivet     PetscMPIInt     rank, size;
4843da7c217SPierre Jolivet 
4853da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
4863da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
4873da7c217SPierre Jolivet     PetscCall(MatGetOwnershipRanges(A, &ranges));
4883da7c217SPierre Jolivet     rStart = ranges[rank];
4893da7c217SPierre Jolivet     rEnd   = ranges[rank + 1];
4909566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
4914325cce7SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
4924325cce7SMatthew G Knepley       PetscInt ncols;
4934325cce7SMatthew G Knepley 
4949566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, r, &ncols, NULL, NULL));
4955f80ce2aSJacob Faibussowitsch       colMax = PetscMax(colMax, ncols);
4969566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL));
4974325cce7SMatthew G Knepley     }
4983da7c217SPierre Jolivet     maxRows = 0;
4993da7c217SPierre Jolivet     for (r = 0; r < size; ++r) maxRows = PetscMax(maxRows, ranges[r + 1] - ranges[r]);
5003da7c217SPierre Jolivet     PetscCall(PetscCalloc2(colMax, &newCols, colMax, &newVals));
5019566063dSJacob Faibussowitsch     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */
5029566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
503cf5d3338SPierre Jolivet     /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd()             */
504cf5d3338SPierre Jolivet     /* that are potentially called many times depending on the distribution of A */
5054325cce7SMatthew G Knepley     for (r = rStart; r < rStart + maxRows; ++r) {
5063da7c217SPierre Jolivet       if (r < rEnd) {
5074325cce7SMatthew G Knepley         const PetscScalar *vals;
5084325cce7SMatthew G Knepley         const PetscInt    *cols;
5093da7c217SPierre Jolivet         PetscInt           ncols, newcols = 0, c;
5104325cce7SMatthew G Knepley 
5119566063dSJacob Faibussowitsch         PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
51234bf4ff8Smarkadams4         nnz0 += ncols - 1;
5134325cce7SMatthew G Knepley         for (c = 0; c < ncols; ++c) {
51482fa6941SPierre Jolivet           if (PetscUnlikely(PetscAbsScalar(vals[c]) <= tol)) newCols[newcols++] = cols[c];
5154325cce7SMatthew G Knepley         }
51634bf4ff8Smarkadams4         nnz1 += ncols - newcols - 1;
5179566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
5189566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
5194325cce7SMatthew G Knepley       }
5209566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5219566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5224325cce7SMatthew G Knepley     }
5239566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
5249566063dSJacob Faibussowitsch     PetscCall(PetscFree2(newCols, newVals));
5259566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */
52634bf4ff8Smarkadams4     if (nnz0 > 0) PetscCall(PetscInfo(NULL, "Filtering left %g%% edges in graph\n", 100 * (double)nnz1 / (double)nnz0));
5275afa97f7SPierre Jolivet     else PetscCall(PetscInfo(NULL, "Warning: %" PetscInt_FMT " edges to filter with %" PetscInt_FMT " rows\n", nnz0, maxRows));
528038df967SPierre Jolivet   }
5292ce66baaSPierre Jolivet   if (compress && A->ops->eliminatezeros) {
5302ce66baaSPierre Jolivet     Mat       B;
53124c7e8a4SPierre Jolivet     PetscBool flg;
5322ce66baaSPierre Jolivet 
53324c7e8a4SPierre Jolivet     PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJHIPSPARSE, MATMPIAIJHIPSPARSE, ""));
53424c7e8a4SPierre Jolivet     if (!flg) {
5352ce66baaSPierre Jolivet       PetscCall(MatEliminateZeros(A, keep));
5362ce66baaSPierre Jolivet       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
5372ce66baaSPierre Jolivet       PetscCall(MatHeaderReplace(A, &B));
5382ce66baaSPierre Jolivet     }
53924c7e8a4SPierre Jolivet   }
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5414325cce7SMatthew G Knepley }
542