xref: /petsc/src/mat/utils/axpy.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
16f79c3a4SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
36f79c3a4SBarry Smith 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTransposeAXPY_Private(Mat Y, PetscScalar a, Mat X, MatStructure str, Mat T)
5d71ae5a4SJacob Faibussowitsch {
6d54c938cSPierre Jolivet   Mat A, F;
75f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat, Mat *);
8d54c938cSPierre Jolivet 
9d54c938cSPierre Jolivet   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)T, "MatTransposeGetMat_C", &f));
11d54c938cSPierre Jolivet   if (f) {
129566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(T, &A));
13d54c938cSPierre Jolivet     if (T == X) {
14013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Explicitly transposing X of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
159566063dSJacob Faibussowitsch       PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &F));
16d54c938cSPierre Jolivet       A = Y;
17d54c938cSPierre Jolivet     } else {
18013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Transposing X because Y of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
199566063dSJacob Faibussowitsch       PetscCall(MatTranspose(X, MAT_INITIAL_MATRIX, &F));
20d54c938cSPierre Jolivet     }
21d54c938cSPierre Jolivet   } else {
229566063dSJacob Faibussowitsch     PetscCall(MatHermitianTransposeGetMat(T, &A));
23d54c938cSPierre Jolivet     if (T == X) {
24013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Explicitly Hermitian transposing X of type MATHERITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
259566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &F));
26d54c938cSPierre Jolivet       A = Y;
27d54c938cSPierre Jolivet     } else {
28013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Hermitian transposing X because Y of type MATHERITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
299566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(X, MAT_INITIAL_MATRIX, &F));
30d54c938cSPierre Jolivet     }
31d54c938cSPierre Jolivet   }
329566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A, a, F, str));
339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35d54c938cSPierre Jolivet }
36d54c938cSPierre Jolivet 
3706be10caSBarry Smith /*@
3821c89e3eSBarry Smith    MatAXPY - Computes Y = a*X + Y.
396f79c3a4SBarry Smith 
40c3339decSBarry Smith    Logically Collective
41fee21e36SBarry Smith 
4298a79cdbSBarry Smith    Input Parameters:
43607cd303SBarry Smith +  a - the scalar multiplier
44607cd303SBarry Smith .  X - the first matrix
45607cd303SBarry Smith .  Y - the second matrix
4627430b45SBarry 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)
47edf9a46cSStefano Zampini 
482860a424SLois Curfman McInnes    Level: intermediate
492860a424SLois Curfman McInnes 
50*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAYPX()`
5106be10caSBarry Smith  @*/
52d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str)
53d71ae5a4SJacob Faibussowitsch {
54646531bbSStefano Zampini   PetscInt  M1, M2, N1, N2;
55c1ac3661SBarry Smith   PetscInt  m1, m2, n1, n2;
56a683ea4aSPierre Jolivet   PetscBool sametype, transpose;
576f79c3a4SBarry Smith 
583a40ed3dSBarry Smith   PetscFunctionBegin;
590700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
60c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y, a, 2);
61646531bbSStefano Zampini   PetscValidHeaderSpecific(X, MAT_CLASSID, 3);
62646531bbSStefano Zampini   PetscCheckSameComm(Y, 1, X, 3);
639566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M1, &N1));
649566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M2, &N2));
659566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m1, &n1));
669566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m2, &n2));
672c71b3e2SJacob 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);
682c71b3e2SJacob 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);
692c71b3e2SJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)");
702c71b3e2SJacob Faibussowitsch   PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)");
713ba16761SJacob Faibussowitsch   if (a == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
72edf9a46cSStefano Zampini   if (Y == X) {
739566063dSJacob Faibussowitsch     PetscCall(MatScale(Y, 1.0 + a));
743ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
75edf9a46cSStefano Zampini   }
76013e2dc7SBarry Smith   PetscCall(PetscObjectObjectTypeCompare((PetscObject)X, (PetscObject)Y, &sametype));
779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0));
7893c87e65SStefano Zampini   if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) {
79dbbe0bcdSBarry Smith     PetscUseTypeMethod(Y, axpy, a, X, str);
80d4bb536fSBarry Smith   } else {
81013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
82a683ea4aSPierre Jolivet     if (transpose) {
839566063dSJacob Faibussowitsch       PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X));
84a683ea4aSPierre Jolivet     } else {
85013e2dc7SBarry Smith       PetscCall(PetscObjectTypeCompareAny((PetscObject)Y, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
86a683ea4aSPierre Jolivet       if (transpose) {
879566063dSJacob Faibussowitsch         PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y));
88a683ea4aSPierre Jolivet       } else {
899566063dSJacob Faibussowitsch         PetscCall(MatAXPY_Basic(Y, a, X, str));
90607cd303SBarry Smith       }
91a683ea4aSPierre Jolivet     }
92a683ea4aSPierre Jolivet   }
939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95607cd303SBarry Smith }
96607cd303SBarry Smith 
97d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
98d71ae5a4SJacob Faibussowitsch {
9927afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL;
100646531bbSStefano Zampini 
101646531bbSStefano Zampini   PetscFunctionBegin;
102646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
1039566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall));
104646531bbSStefano Zampini   if (preall) {
1059566063dSJacob Faibussowitsch     PetscCall((*preall)(Y, X, B));
106646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
107646531bbSStefano Zampini     Mat      preallocator;
108646531bbSStefano Zampini     PetscInt r, rstart, rend;
109646531bbSStefano Zampini     PetscInt m, n, M, N;
110646531bbSStefano Zampini 
1119566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
1129566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
1139566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Y, &M, &N));
1149566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Y, &m, &n));
1159566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator));
1169566063dSJacob Faibussowitsch     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1179566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap));
1189566063dSJacob Faibussowitsch     PetscCall(MatSetUp(preallocator));
1199566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
120646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
121646531bbSStefano Zampini       PetscInt           ncols;
122646531bbSStefano Zampini       const PetscInt    *row;
123646531bbSStefano Zampini       const PetscScalar *vals;
124646531bbSStefano Zampini 
1259566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, r, &ncols, &row, &vals));
1269566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1279566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals));
1289566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, r, &ncols, &row, &vals));
1299566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1309566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals));
131646531bbSStefano Zampini     }
1329566063dSJacob Faibussowitsch     PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1339566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1349566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
1359566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
1369566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
137646531bbSStefano Zampini 
1389566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B));
1399566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name));
1409566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap));
1419566063dSJacob Faibussowitsch     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B));
1429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&preallocator));
143646531bbSStefano Zampini   }
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145646531bbSStefano Zampini }
146646531bbSStefano Zampini 
147d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str)
148d71ae5a4SJacob Faibussowitsch {
149be705e3aSPierre Jolivet   PetscBool isshell, isdense, isnest;
150edf9a46cSStefano Zampini 
151edf9a46cSStefano Zampini   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(MatIsShell(Y, &isshell));
153edf9a46cSStefano Zampini   if (isshell) { /* MatShell has special support for AXPY */
154edf9a46cSStefano Zampini     PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure);
155edf9a46cSStefano Zampini 
1569566063dSJacob Faibussowitsch     PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f));
157edf9a46cSStefano Zampini     if (f) {
1589566063dSJacob Faibussowitsch       PetscCall((*f)(Y, a, X, str));
1593ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
160edf9a46cSStefano Zampini     }
161edf9a46cSStefano Zampini   }
16293c87e65SStefano Zampini   /* no need to preallocate if Y is dense */
1639566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, ""));
164be705e3aSPierre Jolivet   if (isdense) {
1659566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &isnest));
166be705e3aSPierre Jolivet     if (isnest) {
1679566063dSJacob Faibussowitsch       PetscCall(MatAXPY_Dense_Nest(Y, a, X));
1683ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
169be705e3aSPierre Jolivet     }
170be705e3aSPierre Jolivet     if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN;
171be705e3aSPierre Jolivet   }
172a6ab3590SBarry Smith   if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) {
173edf9a46cSStefano Zampini     PetscInt           i, start, end, j, ncols, m, n;
17438baddfdSBarry Smith     const PetscInt    *row;
175b3cc6726SBarry Smith     PetscScalar       *val;
176b3cc6726SBarry Smith     const PetscScalar *vals;
177607cd303SBarry Smith 
1789566063dSJacob Faibussowitsch     PetscCall(MatGetSize(X, &m, &n));
1799566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(X, &start, &end));
1809566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
181f4df32b1SMatthew Knepley     if (a == 1.0) {
182d4bb536fSBarry Smith       for (i = start; i < end; i++) {
1839566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
1849566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES));
1859566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
186d4bb536fSBarry Smith       }
187d4bb536fSBarry Smith     } else {
18821a3365eSStefano Zampini       PetscInt vs = 100;
18921a3365eSStefano Zampini       /* realloc if needed, as this function may be used in parallel */
1909566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(vs, &val));
19106be10caSBarry Smith       for (i = start; i < end; i++) {
1929566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
19321a3365eSStefano Zampini         if (vs < ncols) {
19421a3365eSStefano Zampini           vs = PetscMin(2 * ncols, n);
1959566063dSJacob Faibussowitsch           PetscCall(PetscRealloc(vs * sizeof(*val), &val));
1966f79c3a4SBarry Smith         }
19721a3365eSStefano Zampini         for (j = 0; j < ncols; j++) val[j] = a * vals[j];
1989566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES));
1999566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
2006f79c3a4SBarry Smith       }
2019566063dSJacob Faibussowitsch       PetscCall(PetscFree(val));
202d4bb536fSBarry Smith     }
2039566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
2049566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
2059566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
206edf9a46cSStefano Zampini   } else {
207edf9a46cSStefano Zampini     Mat B;
208edf9a46cSStefano Zampini 
2099566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B));
2109566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
2119566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y, &B));
212edf9a46cSStefano Zampini   }
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2146f79c3a4SBarry Smith }
215052efed2SBarry Smith 
216d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str)
217d71ae5a4SJacob Faibussowitsch {
218ec7775f6SShri Abhyankar   PetscInt           i, start, end, j, ncols, m, n;
219ec7775f6SShri Abhyankar   const PetscInt    *row;
220ec7775f6SShri Abhyankar   PetscScalar       *val;
221ec7775f6SShri Abhyankar   const PetscScalar *vals;
222ec7775f6SShri Abhyankar 
223ec7775f6SShri Abhyankar   PetscFunctionBegin;
2249566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &m, &n));
2259566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(X, &start, &end));
2269566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(Y));
2279566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(X));
228ec7775f6SShri Abhyankar   if (a == 1.0) {
229ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2309566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2319566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2329566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
233ec7775f6SShri Abhyankar 
2349566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2359566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2369566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
237ec7775f6SShri Abhyankar     }
238ec7775f6SShri Abhyankar   } else {
23921a3365eSStefano Zampini     PetscInt vs = 100;
24021a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
2419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(vs, &val));
242ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2439566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2449566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2459566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
246ec7775f6SShri Abhyankar 
2479566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
24821a3365eSStefano Zampini       if (vs < ncols) {
24921a3365eSStefano Zampini         vs = PetscMin(2 * ncols, n);
2509566063dSJacob Faibussowitsch         PetscCall(PetscRealloc(vs * sizeof(*val), &val));
251ec7775f6SShri Abhyankar       }
25221a3365eSStefano Zampini       for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2539566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES));
2549566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
255ec7775f6SShri Abhyankar     }
2569566063dSJacob Faibussowitsch     PetscCall(PetscFree(val));
257ec7775f6SShri Abhyankar   }
2589566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(Y));
2599566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(X));
2609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263ec7775f6SShri Abhyankar }
264ec7775f6SShri Abhyankar 
265052efed2SBarry Smith /*@
2662ef1f0ffSBarry Smith    MatShift - Computes `Y =  Y + a I`, where `a` is a `PetscScalar`
267052efed2SBarry Smith 
268c3339decSBarry Smith    Neighbor-wise Collective
269fee21e36SBarry Smith 
27098a79cdbSBarry Smith    Input Parameters:
27127430b45SBarry Smith +  Y - the matrix
27227430b45SBarry Smith -  a - the `PetscScalar`
27398a79cdbSBarry Smith 
2742860a424SLois Curfman McInnes    Level: intermediate
2752860a424SLois Curfman McInnes 
27695452b02SPatrick Sanan    Notes:
27727430b45SBarry 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)
27884e19e3eSJunchao Zhang 
2792ef1f0ffSBarry Smith     If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2806f33a894SBarry 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
281edf9a46cSStefano Zampini    entry). No operation is performed when a is zero.
2826f33a894SBarry Smith 
28327430b45SBarry Smith    To form Y = Y + diag(V) use `MatDiagonalSet()`
2840c0fd78eSBarry Smith 
285*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()`
286052efed2SBarry Smith  @*/
287d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift(Mat Y, PetscScalar a)
288d71ae5a4SJacob Faibussowitsch {
2893a40ed3dSBarry Smith   PetscFunctionBegin;
2900700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
2915f80ce2aSJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2925f80ce2aSJacob Faibussowitsch   PetscCheck(!Y->factortype, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2934994cf47SJed Brown   MatCheckPreallocated(Y, 1);
2943ba16761SJacob Faibussowitsch   if (a == 0.0) PetscFunctionReturn(PETSC_SUCCESS);
295b50b34bdSBarry Smith 
296dbbe0bcdSBarry Smith   if (Y->ops->shift) PetscUseTypeMethod(Y, shift, a);
2979566063dSJacob Faibussowitsch   else PetscCall(MatShift_Basic(Y, a));
2987d68702bSBarry Smith 
2999566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301052efed2SBarry Smith }
3026d84be18SBarry Smith 
303d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is)
304d71ae5a4SJacob Faibussowitsch {
30567576b8bSBarry Smith   PetscInt           i, start, end;
306fa2e578bSStefano Zampini   const PetscScalar *v;
30709f38230SBarry Smith 
30809f38230SBarry Smith   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y, &start, &end));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(D, &v));
31148a46eb9SPierre Jolivet   for (i = start; i < end; i++) PetscCall(MatSetValues(Y, 1, &i, 1, &i, v + i - start, is));
3129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(D, &v));
3139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
3149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31609f38230SBarry Smith }
31709f38230SBarry Smith 
3186d84be18SBarry Smith /*@
3192ef1f0ffSBarry Smith    MatDiagonalSet - Computes `Y` = `Y` + `D`, where `D` is a diagonal matrix
3202ef1f0ffSBarry Smith    that is represented as a vector. Or Y[i,i] = D[i] if `InsertMode` is
3212ef1f0ffSBarry Smith    `INSERT_VALUES`.
3226d84be18SBarry Smith 
323c3339decSBarry Smith    Neighbor-wise Collective
32448e586daSJose E. Roman 
3256d84be18SBarry Smith    Input Parameters:
32698a79cdbSBarry Smith +  Y - the input matrix
327f56f2b3fSBarry Smith .  D - the diagonal matrix, represented as a vector
3282ef1f0ffSBarry Smith -  i - `INSERT_VALUES` or `ADD_VALUES`
3296f33a894SBarry Smith 
3302860a424SLois Curfman McInnes    Level: intermediate
3312860a424SLois Curfman McInnes 
3322ef1f0ffSBarry Smith    Note:
3332ef1f0ffSBarry Smith     If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3342ef1f0ffSBarry 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
3352ef1f0ffSBarry Smith    entry).
3362ef1f0ffSBarry Smith 
337*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatShift()`, `MatScale()`, `MatDiagonalScale()`
3386d84be18SBarry Smith @*/
339d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is)
340d71ae5a4SJacob Faibussowitsch {
34167576b8bSBarry Smith   PetscInt matlocal, veclocal;
3426d84be18SBarry Smith 
3433a40ed3dSBarry Smith   PetscFunctionBegin;
3440700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3450700a824SBarry Smith   PetscValidHeaderSpecific(D, VEC_CLASSID, 2);
3469566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &matlocal, NULL));
3479566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(D, &veclocal));
3485f80ce2aSJacob 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);
349dbbe0bcdSBarry Smith   if (Y->ops->diagonalset) PetscUseTypeMethod(Y, diagonalset, D, is);
350dbbe0bcdSBarry Smith   else PetscCall(MatDiagonalSet_Default(Y, D, is));
3519566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3536d84be18SBarry Smith }
354d4bb536fSBarry Smith 
355d4bb536fSBarry Smith /*@
35604aac2b0SHong Zhang    MatAYPX - Computes Y = a*Y + X.
357d4bb536fSBarry Smith 
3582ef1f0ffSBarry Smith    Logically Collective
359fee21e36SBarry Smith 
36098a79cdbSBarry Smith    Input Parameters:
3612ef1f0ffSBarry Smith +  a - the `PetscScalar` multiplier
36204aac2b0SHong Zhang .  Y - the first matrix
36304aac2b0SHong Zhang .  X - the second matrix
364aaa8cc7dSPierre 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)
365d4bb536fSBarry Smith 
3662860a424SLois Curfman McInnes    Level: intermediate
3672860a424SLois Curfman McInnes 
368*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAXPY()`
369d4bb536fSBarry Smith  @*/
370d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAYPX(Mat Y, PetscScalar a, Mat X, MatStructure str)
371d71ae5a4SJacob Faibussowitsch {
3723a40ed3dSBarry Smith   PetscFunctionBegin;
3739566063dSJacob Faibussowitsch   PetscCall(MatScale(Y, a));
3749566063dSJacob Faibussowitsch   PetscCall(MatAXPY(Y, 1.0, X, str));
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
376d4bb536fSBarry Smith }
377b0a32e0cSBarry Smith 
378b0a32e0cSBarry Smith /*@
3790bacdadaSStefano Zampini     MatComputeOperator - Computes the explicit matrix
380b0a32e0cSBarry Smith 
381c3339decSBarry Smith     Collective
382b0a32e0cSBarry Smith 
383d8d19677SJose E. Roman     Input Parameters:
384186a3e20SStefano Zampini +   inmat - the matrix
385186a3e20SStefano Zampini -   mattype - the matrix type for the explicit operator
386b0a32e0cSBarry Smith 
387b0a32e0cSBarry Smith     Output Parameter:
388a5b23f4aSJose E. Roman .   mat - the explicit  operator
389b0a32e0cSBarry Smith 
3902ef1f0ffSBarry Smith     Level: advanced
3912ef1f0ffSBarry Smith 
39211a5261eSBarry Smith     Note:
393186a3e20SStefano Zampini     This computation is done by applying the operators to columns of the identity matrix.
394186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
3952ef1f0ffSBarry Smith     Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
396b0a32e0cSBarry Smith 
397*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperatorTranspose()`
398b0a32e0cSBarry Smith @*/
399d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperator(Mat inmat, MatType mattype, Mat *mat)
400d71ae5a4SJacob Faibussowitsch {
401b0a32e0cSBarry Smith   PetscFunctionBegin;
4020700a824SBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
403186a3e20SStefano Zampini   PetscValidPointer(mat, 3);
4049566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(inmat, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40624f910e3SHong Zhang }
4074325cce7SMatthew G Knepley 
4084325cce7SMatthew G Knepley /*@
4090bacdadaSStefano Zampini     MatComputeOperatorTranspose - Computes the explicit matrix representation of
4102ef1f0ffSBarry Smith         a give matrix that can apply `MatMultTranspose()`
411f3b1f45cSBarry Smith 
412c3339decSBarry Smith     Collective
413f3b1f45cSBarry Smith 
4146b867d5aSJose E. Roman     Input Parameters:
4156b867d5aSJose E. Roman +   inmat - the matrix
4166b867d5aSJose E. Roman -   mattype - the matrix type for the explicit operator
417f3b1f45cSBarry Smith 
418f3b1f45cSBarry Smith     Output Parameter:
419a5b23f4aSJose E. Roman .   mat - the explicit  operator transposed
420f3b1f45cSBarry Smith 
4212ef1f0ffSBarry Smith     Level: advanced
4222ef1f0ffSBarry Smith 
42311a5261eSBarry Smith     Note:
424186a3e20SStefano Zampini     This computation is done by applying the transpose of the operator to columns of the identity matrix.
425186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
4262ef1f0ffSBarry Smith     Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
427f3b1f45cSBarry Smith 
428*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperator()`
429f3b1f45cSBarry Smith @*/
430d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperatorTranspose(Mat inmat, MatType mattype, Mat *mat)
431d71ae5a4SJacob Faibussowitsch {
432186a3e20SStefano Zampini   Mat A;
433f3b1f45cSBarry Smith 
434f3b1f45cSBarry Smith   PetscFunctionBegin;
435f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
436186a3e20SStefano Zampini   PetscValidPointer(mat, 3);
4379566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(inmat, &A));
4389566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(A, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
441f3b1f45cSBarry Smith }
442f3b1f45cSBarry Smith 
443f3b1f45cSBarry Smith /*@
4444325cce7SMatthew G Knepley   MatChop - Set all values in the matrix less than the tolerance to zero
4454325cce7SMatthew G Knepley 
4464325cce7SMatthew G Knepley   Input Parameters:
4474325cce7SMatthew G Knepley + A   - The matrix
4484325cce7SMatthew G Knepley - tol - The zero tolerance
4494325cce7SMatthew G Knepley 
4504325cce7SMatthew G Knepley   Level: intermediate
4514325cce7SMatthew G Knepley 
452*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatZeroEntries()`
4533fc99919SSatish Balay  @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode MatChop(Mat A, PetscReal tol)
455d71ae5a4SJacob Faibussowitsch {
456038df967SPierre Jolivet   Mat          a;
4574325cce7SMatthew G Knepley   PetscScalar *newVals;
458cf5d3338SPierre Jolivet   PetscInt    *newCols, rStart, rEnd, numRows, maxRows, r, colMax = 0;
459038df967SPierre Jolivet   PetscBool    flg;
4604325cce7SMatthew G Knepley 
4614325cce7SMatthew G Knepley   PetscFunctionBegin;
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, ""));
463038df967SPierre Jolivet   if (flg) {
4649566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLocalMatrix(A, &a));
4659566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(a, &r));
4669566063dSJacob Faibussowitsch     PetscCall(MatGetSize(a, &rStart, &rEnd));
4679566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(a, &newVals));
468038df967SPierre Jolivet     for (; colMax < rEnd; ++colMax) {
469ad540459SPierre Jolivet       for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) < tol ? 0.0 : newVals[maxRows + colMax * r];
470038df967SPierre Jolivet     }
4719566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(a, &newVals));
472038df967SPierre Jolivet   } else {
4739566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
4749566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
4754325cce7SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
4764325cce7SMatthew G Knepley       PetscInt ncols;
4774325cce7SMatthew G Knepley 
4789566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, r, &ncols, NULL, NULL));
4795f80ce2aSJacob Faibussowitsch       colMax = PetscMax(colMax, ncols);
4809566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL));
4814325cce7SMatthew G Knepley     }
4824325cce7SMatthew G Knepley     numRows = rEnd - rStart;
4831c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A)));
4849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(colMax, &newCols, colMax, &newVals));
4859566063dSJacob Faibussowitsch     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */
4869566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
487cf5d3338SPierre Jolivet     /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd()             */
488cf5d3338SPierre Jolivet     /* that are potentially called many times depending on the distribution of A */
4894325cce7SMatthew G Knepley     for (r = rStart; r < rStart + maxRows; ++r) {
4904325cce7SMatthew G Knepley       const PetscScalar *vals;
4914325cce7SMatthew G Knepley       const PetscInt    *cols;
492fad45679SMatthew G. Knepley       PetscInt           ncols, newcols, c;
4934325cce7SMatthew G Knepley 
4944325cce7SMatthew G Knepley       if (r < rEnd) {
4959566063dSJacob Faibussowitsch         PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
4964325cce7SMatthew G Knepley         for (c = 0; c < ncols; ++c) {
4974325cce7SMatthew G Knepley           newCols[c] = cols[c];
4984325cce7SMatthew G Knepley           newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
4994325cce7SMatthew G Knepley         }
500fad45679SMatthew G. Knepley         newcols = ncols;
5019566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
5029566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
5034325cce7SMatthew G Knepley       }
5049566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5059566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5064325cce7SMatthew G Knepley     }
5079566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
5089566063dSJacob Faibussowitsch     PetscCall(PetscFree2(newCols, newVals));
5099566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */
510038df967SPierre Jolivet   }
5113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5124325cce7SMatthew G Knepley }
513