xref: /petsc/src/mat/utils/axpy.c (revision f13dfd9ea68e0ddeee984e65c377a1819eab8a8a)
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;
65f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat, Mat *);
7d54c938cSPierre Jolivet 
8d54c938cSPierre Jolivet   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)T, "MatTransposeGetMat_C", &f));
10d54c938cSPierre Jolivet   if (f) {
119566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(T, &A));
12d54c938cSPierre Jolivet     if (T == X) {
13013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Explicitly transposing X of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
149566063dSJacob Faibussowitsch       PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &F));
15d54c938cSPierre Jolivet       A = Y;
16d54c938cSPierre Jolivet     } else {
17013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Transposing X because Y of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
189566063dSJacob Faibussowitsch       PetscCall(MatTranspose(X, MAT_INITIAL_MATRIX, &F));
19d54c938cSPierre Jolivet     }
20d54c938cSPierre Jolivet   } else {
219566063dSJacob Faibussowitsch     PetscCall(MatHermitianTransposeGetMat(T, &A));
22d54c938cSPierre Jolivet     if (T == X) {
2352cd20c4SPierre Jolivet       PetscCall(PetscInfo(NULL, "Explicitly Hermitian transposing X of type MATHERMITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
249566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &F));
25d54c938cSPierre Jolivet       A = Y;
26d54c938cSPierre Jolivet     } else {
2752cd20c4SPierre Jolivet       PetscCall(PetscInfo(NULL, "Hermitian transposing X because Y of type MATHERMITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
289566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(X, MAT_INITIAL_MATRIX, &F));
29d54c938cSPierre Jolivet     }
30d54c938cSPierre Jolivet   }
319566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A, a, F, str));
329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34d54c938cSPierre Jolivet }
35d54c938cSPierre Jolivet 
3606be10caSBarry Smith /*@
3721c89e3eSBarry Smith   MatAXPY - Computes Y = a*X + Y.
386f79c3a4SBarry Smith 
39c3339decSBarry Smith   Logically Collective
40fee21e36SBarry Smith 
4198a79cdbSBarry Smith   Input Parameters:
42607cd303SBarry Smith + a   - the scalar multiplier
43607cd303SBarry Smith . X   - the first matrix
44607cd303SBarry Smith . Y   - the second matrix
4527430b45SBarry 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)
46edf9a46cSStefano Zampini 
472860a424SLois Curfman McInnes   Level: intermediate
482860a424SLois Curfman McInnes 
491cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAYPX()`
5006be10caSBarry Smith  @*/
51d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str)
52d71ae5a4SJacob Faibussowitsch {
53646531bbSStefano Zampini   PetscInt  M1, M2, N1, N2;
54c1ac3661SBarry Smith   PetscInt  m1, m2, n1, n2;
55a683ea4aSPierre Jolivet   PetscBool sametype, transpose;
566f79c3a4SBarry Smith 
573a40ed3dSBarry Smith   PetscFunctionBegin;
580700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
59c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y, a, 2);
60646531bbSStefano Zampini   PetscValidHeaderSpecific(X, MAT_CLASSID, 3);
61646531bbSStefano Zampini   PetscCheckSameComm(Y, 1, X, 3);
629566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M1, &N1));
639566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M2, &N2));
649566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m1, &n1));
659566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m2, &n2));
662c71b3e2SJacob 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);
672c71b3e2SJacob 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);
682c71b3e2SJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)");
692c71b3e2SJacob Faibussowitsch   PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)");
703ba16761SJacob Faibussowitsch   if (a == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
71edf9a46cSStefano Zampini   if (Y == X) {
729566063dSJacob Faibussowitsch     PetscCall(MatScale(Y, 1.0 + a));
733ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
74edf9a46cSStefano Zampini   }
75013e2dc7SBarry Smith   PetscCall(PetscObjectObjectTypeCompare((PetscObject)X, (PetscObject)Y, &sametype));
769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0));
7793c87e65SStefano Zampini   if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) {
78dbbe0bcdSBarry Smith     PetscUseTypeMethod(Y, axpy, a, X, str);
79d4bb536fSBarry Smith   } else {
80013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
81a683ea4aSPierre Jolivet     if (transpose) {
829566063dSJacob Faibussowitsch       PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X));
83a683ea4aSPierre Jolivet     } else {
84013e2dc7SBarry Smith       PetscCall(PetscObjectTypeCompareAny((PetscObject)Y, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
85a683ea4aSPierre Jolivet       if (transpose) {
869566063dSJacob Faibussowitsch         PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y));
87a683ea4aSPierre Jolivet       } else {
889566063dSJacob Faibussowitsch         PetscCall(MatAXPY_Basic(Y, a, X, str));
89607cd303SBarry Smith       }
90a683ea4aSPierre Jolivet     }
91a683ea4aSPierre Jolivet   }
929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0));
933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
94607cd303SBarry Smith }
95607cd303SBarry Smith 
96d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
97d71ae5a4SJacob Faibussowitsch {
9827afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL;
99646531bbSStefano Zampini 
100646531bbSStefano Zampini   PetscFunctionBegin;
101646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall));
103646531bbSStefano Zampini   if (preall) {
1049566063dSJacob Faibussowitsch     PetscCall((*preall)(Y, X, B));
105646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
106646531bbSStefano Zampini     Mat      preallocator;
107646531bbSStefano Zampini     PetscInt r, rstart, rend;
108646531bbSStefano Zampini     PetscInt m, n, M, N;
109646531bbSStefano Zampini 
1109566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
1119566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
1129566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Y, &M, &N));
1139566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Y, &m, &n));
1149566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator));
1159566063dSJacob Faibussowitsch     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1169566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap));
1179566063dSJacob Faibussowitsch     PetscCall(MatSetUp(preallocator));
1189566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
119646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
120646531bbSStefano Zampini       PetscInt           ncols;
121646531bbSStefano Zampini       const PetscInt    *row;
122646531bbSStefano Zampini       const PetscScalar *vals;
123646531bbSStefano Zampini 
1249566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, r, &ncols, &row, &vals));
1259566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1269566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals));
1279566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, r, &ncols, &row, &vals));
1289566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1299566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals));
130646531bbSStefano Zampini     }
1319566063dSJacob Faibussowitsch     PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1339566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
1349566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
1359566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
136646531bbSStefano Zampini 
1379566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B));
1389566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name));
1399566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap));
1409566063dSJacob Faibussowitsch     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B));
1419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&preallocator));
142646531bbSStefano Zampini   }
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144646531bbSStefano Zampini }
145646531bbSStefano Zampini 
146d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str)
147d71ae5a4SJacob Faibussowitsch {
148be705e3aSPierre Jolivet   PetscBool isshell, isdense, isnest;
149edf9a46cSStefano Zampini 
150edf9a46cSStefano Zampini   PetscFunctionBegin;
1519566063dSJacob Faibussowitsch   PetscCall(MatIsShell(Y, &isshell));
152edf9a46cSStefano Zampini   if (isshell) { /* MatShell has special support for AXPY */
153edf9a46cSStefano Zampini     PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure);
154edf9a46cSStefano Zampini 
1559566063dSJacob Faibussowitsch     PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f));
156edf9a46cSStefano Zampini     if (f) {
1579566063dSJacob Faibussowitsch       PetscCall((*f)(Y, a, X, str));
1583ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
159edf9a46cSStefano Zampini     }
160edf9a46cSStefano Zampini   }
16193c87e65SStefano Zampini   /* no need to preallocate if Y is dense */
1629566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, ""));
163be705e3aSPierre Jolivet   if (isdense) {
1649566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &isnest));
165be705e3aSPierre Jolivet     if (isnest) {
1669566063dSJacob Faibussowitsch       PetscCall(MatAXPY_Dense_Nest(Y, a, X));
1673ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
168be705e3aSPierre Jolivet     }
169be705e3aSPierre Jolivet     if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN;
170be705e3aSPierre Jolivet   }
171a6ab3590SBarry Smith   if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) {
172edf9a46cSStefano Zampini     PetscInt           i, start, end, j, ncols, m, n;
17338baddfdSBarry Smith     const PetscInt    *row;
174b3cc6726SBarry Smith     PetscScalar       *val;
175b3cc6726SBarry Smith     const PetscScalar *vals;
176607cd303SBarry Smith 
1779566063dSJacob Faibussowitsch     PetscCall(MatGetSize(X, &m, &n));
1789566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(X, &start, &end));
1799566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
180f4df32b1SMatthew Knepley     if (a == 1.0) {
181d4bb536fSBarry Smith       for (i = start; i < end; i++) {
1829566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
1839566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES));
1849566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
185d4bb536fSBarry Smith       }
186d4bb536fSBarry Smith     } else {
18721a3365eSStefano Zampini       PetscInt vs = 100;
18821a3365eSStefano Zampini       /* realloc if needed, as this function may be used in parallel */
1899566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(vs, &val));
19006be10caSBarry Smith       for (i = start; i < end; i++) {
1919566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
19221a3365eSStefano Zampini         if (vs < ncols) {
19321a3365eSStefano Zampini           vs = PetscMin(2 * ncols, n);
1949566063dSJacob Faibussowitsch           PetscCall(PetscRealloc(vs * sizeof(*val), &val));
1956f79c3a4SBarry Smith         }
19621a3365eSStefano Zampini         for (j = 0; j < ncols; j++) val[j] = a * vals[j];
1979566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES));
1989566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
1996f79c3a4SBarry Smith       }
2009566063dSJacob Faibussowitsch       PetscCall(PetscFree(val));
201d4bb536fSBarry Smith     }
2029566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
2039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
2049566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
205edf9a46cSStefano Zampini   } else {
206edf9a46cSStefano Zampini     Mat B;
207edf9a46cSStefano Zampini 
2089566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B));
2099566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
2109566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y, &B));
211edf9a46cSStefano Zampini   }
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2136f79c3a4SBarry Smith }
214052efed2SBarry Smith 
215d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str)
216d71ae5a4SJacob Faibussowitsch {
217ec7775f6SShri Abhyankar   PetscInt           i, start, end, j, ncols, m, n;
218ec7775f6SShri Abhyankar   const PetscInt    *row;
219ec7775f6SShri Abhyankar   PetscScalar       *val;
220ec7775f6SShri Abhyankar   const PetscScalar *vals;
221ec7775f6SShri Abhyankar 
222ec7775f6SShri Abhyankar   PetscFunctionBegin;
2239566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &m, &n));
2249566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(X, &start, &end));
2259566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(Y));
2269566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(X));
227ec7775f6SShri Abhyankar   if (a == 1.0) {
228ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2299566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2309566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2319566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
232ec7775f6SShri Abhyankar 
2339566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2349566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2359566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
236ec7775f6SShri Abhyankar     }
237ec7775f6SShri Abhyankar   } else {
23821a3365eSStefano Zampini     PetscInt vs = 100;
23921a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
2409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(vs, &val));
241ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2429566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2439566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2449566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
245ec7775f6SShri Abhyankar 
2469566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
24721a3365eSStefano Zampini       if (vs < ncols) {
24821a3365eSStefano Zampini         vs = PetscMin(2 * ncols, n);
2499566063dSJacob Faibussowitsch         PetscCall(PetscRealloc(vs * sizeof(*val), &val));
250ec7775f6SShri Abhyankar       }
25121a3365eSStefano Zampini       for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2529566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES));
2539566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
254ec7775f6SShri Abhyankar     }
2559566063dSJacob Faibussowitsch     PetscCall(PetscFree(val));
256ec7775f6SShri Abhyankar   }
2579566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(Y));
2589566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(X));
2599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262ec7775f6SShri Abhyankar }
263ec7775f6SShri Abhyankar 
264052efed2SBarry Smith /*@
2652ef1f0ffSBarry Smith   MatShift - Computes `Y =  Y + a I`, where `a` is a `PetscScalar`
266052efed2SBarry Smith 
267c3339decSBarry Smith   Neighbor-wise Collective
268fee21e36SBarry Smith 
26998a79cdbSBarry Smith   Input Parameters:
27027430b45SBarry Smith + Y - the matrix
27127430b45SBarry Smith - a - the `PetscScalar`
27298a79cdbSBarry Smith 
2732860a424SLois Curfman McInnes   Level: intermediate
2742860a424SLois Curfman McInnes 
27595452b02SPatrick Sanan   Notes:
27627430b45SBarry 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)
27784e19e3eSJunchao Zhang 
2782ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2796f33a894SBarry 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
280edf9a46cSStefano Zampini   entry). No operation is performed when a is zero.
2816f33a894SBarry Smith 
28227430b45SBarry Smith   To form Y = Y + diag(V) use `MatDiagonalSet()`
2830c0fd78eSBarry Smith 
2841cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()`
285052efed2SBarry Smith  @*/
286d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift(Mat Y, PetscScalar a)
287d71ae5a4SJacob Faibussowitsch {
2883a40ed3dSBarry Smith   PetscFunctionBegin;
2890700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
2905f80ce2aSJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2915f80ce2aSJacob Faibussowitsch   PetscCheck(!Y->factortype, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2924994cf47SJed Brown   MatCheckPreallocated(Y, 1);
2933ba16761SJacob Faibussowitsch   if (a == 0.0) PetscFunctionReturn(PETSC_SUCCESS);
294b50b34bdSBarry Smith 
295dbbe0bcdSBarry Smith   if (Y->ops->shift) PetscUseTypeMethod(Y, shift, a);
2969566063dSJacob Faibussowitsch   else PetscCall(MatShift_Basic(Y, a));
2977d68702bSBarry Smith 
2989566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
300052efed2SBarry Smith }
3016d84be18SBarry Smith 
302d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is)
303d71ae5a4SJacob Faibussowitsch {
30467576b8bSBarry Smith   PetscInt           i, start, end;
305fa2e578bSStefano Zampini   const PetscScalar *v;
30609f38230SBarry Smith 
30709f38230SBarry Smith   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y, &start, &end));
3099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(D, &v));
31048a46eb9SPierre Jolivet   for (i = start; i < end; i++) PetscCall(MatSetValues(Y, 1, &i, 1, &i, v + i - start, is));
3119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(D, &v));
3129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
3139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31509f38230SBarry Smith }
31609f38230SBarry Smith 
3176d84be18SBarry Smith /*@
3182ef1f0ffSBarry Smith   MatDiagonalSet - Computes `Y` = `Y` + `D`, where `D` is a diagonal matrix
3192ef1f0ffSBarry Smith   that is represented as a vector. Or Y[i,i] = D[i] if `InsertMode` is
3202ef1f0ffSBarry Smith   `INSERT_VALUES`.
3216d84be18SBarry Smith 
322c3339decSBarry Smith   Neighbor-wise Collective
32348e586daSJose E. Roman 
3246d84be18SBarry Smith   Input Parameters:
32598a79cdbSBarry Smith + Y  - the input matrix
326f56f2b3fSBarry Smith . D  - the diagonal matrix, represented as a vector
327fe59aa6dSJacob Faibussowitsch - is - `INSERT_VALUES` or `ADD_VALUES`
3286f33a894SBarry Smith 
3292860a424SLois Curfman McInnes   Level: intermediate
3302860a424SLois Curfman McInnes 
3312ef1f0ffSBarry Smith   Note:
3322ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3332ef1f0ffSBarry 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
3342ef1f0ffSBarry Smith   entry).
3352ef1f0ffSBarry Smith 
3361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatShift()`, `MatScale()`, `MatDiagonalScale()`
3376d84be18SBarry Smith @*/
338d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is)
339d71ae5a4SJacob Faibussowitsch {
34067576b8bSBarry Smith   PetscInt matlocal, veclocal;
3416d84be18SBarry Smith 
3423a40ed3dSBarry Smith   PetscFunctionBegin;
3430700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3440700a824SBarry Smith   PetscValidHeaderSpecific(D, VEC_CLASSID, 2);
3452fbb122aSJose E. Roman   MatCheckPreallocated(Y, 1);
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 
3681cc06b55SBarry 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 
378fe59aa6dSJacob Faibussowitsch /*@C
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:
393*f13dfd9eSBarry Smith   This computation is done by applying the operator 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 
3971cc06b55SBarry 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);
4034f572ea9SToby Isaac   PetscAssertPointer(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 
408fe59aa6dSJacob Faibussowitsch /*@C
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 
4281cc06b55SBarry 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);
4364f572ea9SToby Isaac   PetscAssertPointer(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 /*@
4442ce66baaSPierre 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
4454325cce7SMatthew G Knepley 
4464325cce7SMatthew G Knepley   Input Parameters:
4474325cce7SMatthew G Knepley + A        - The matrix
4482ce66baaSPierre Jolivet . tol      - The zero tolerance
4492ce66baaSPierre 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
4502ce66baaSPierre 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
4514325cce7SMatthew G Knepley 
4524325cce7SMatthew G Knepley   Level: intermediate
4534325cce7SMatthew G Knepley 
454d9b70fcdSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatZeroEntries()`, `MatEliminateZeros()`, `VecFilter()`
4553fc99919SSatish Balay  @*/
4562ce66baaSPierre Jolivet PetscErrorCode MatFilter(Mat A, PetscReal tol, PetscBool compress, PetscBool keep)
457d71ae5a4SJacob Faibussowitsch {
458038df967SPierre Jolivet   Mat          a;
4594325cce7SMatthew G Knepley   PetscScalar *newVals;
46034bf4ff8Smarkadams4   PetscInt    *newCols, rStart, rEnd, maxRows, r, colMax = 0, nnz0 = 0, nnz1 = 0;
461038df967SPierre Jolivet   PetscBool    flg;
4624325cce7SMatthew G Knepley 
4634325cce7SMatthew G Knepley   PetscFunctionBegin;
4649566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, ""));
465038df967SPierre Jolivet   if (flg) {
4669566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLocalMatrix(A, &a));
4679566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(a, &r));
4689566063dSJacob Faibussowitsch     PetscCall(MatGetSize(a, &rStart, &rEnd));
4699566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(a, &newVals));
470038df967SPierre Jolivet     for (; colMax < rEnd; ++colMax) {
47182fa6941SPierre Jolivet       for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) <= tol ? 0.0 : newVals[maxRows + colMax * r];
472038df967SPierre Jolivet     }
4739566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(a, &newVals));
474038df967SPierre Jolivet   } else {
4753da7c217SPierre Jolivet     const PetscInt *ranges;
4763da7c217SPierre Jolivet     PetscMPIInt     rank, size;
4773da7c217SPierre Jolivet 
4783da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
4793da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
4803da7c217SPierre Jolivet     PetscCall(MatGetOwnershipRanges(A, &ranges));
4813da7c217SPierre Jolivet     rStart = ranges[rank];
4823da7c217SPierre Jolivet     rEnd   = ranges[rank + 1];
4839566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
4844325cce7SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
4854325cce7SMatthew G Knepley       PetscInt ncols;
4864325cce7SMatthew G Knepley 
4879566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, r, &ncols, NULL, NULL));
4885f80ce2aSJacob Faibussowitsch       colMax = PetscMax(colMax, ncols);
4899566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL));
4904325cce7SMatthew G Knepley     }
4913da7c217SPierre Jolivet     maxRows = 0;
4923da7c217SPierre Jolivet     for (r = 0; r < size; ++r) maxRows = PetscMax(maxRows, ranges[r + 1] - ranges[r]);
4933da7c217SPierre Jolivet     PetscCall(PetscCalloc2(colMax, &newCols, colMax, &newVals));
4949566063dSJacob Faibussowitsch     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */
4959566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
496cf5d3338SPierre Jolivet     /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd()             */
497cf5d3338SPierre Jolivet     /* that are potentially called many times depending on the distribution of A */
4984325cce7SMatthew G Knepley     for (r = rStart; r < rStart + maxRows; ++r) {
4993da7c217SPierre Jolivet       if (r < rEnd) {
5004325cce7SMatthew G Knepley         const PetscScalar *vals;
5014325cce7SMatthew G Knepley         const PetscInt    *cols;
5023da7c217SPierre Jolivet         PetscInt           ncols, newcols = 0, c;
5034325cce7SMatthew G Knepley 
5049566063dSJacob Faibussowitsch         PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
50534bf4ff8Smarkadams4         nnz0 += ncols - 1;
5064325cce7SMatthew G Knepley         for (c = 0; c < ncols; ++c) {
50782fa6941SPierre Jolivet           if (PetscUnlikely(PetscAbsScalar(vals[c]) <= tol)) newCols[newcols++] = cols[c];
5084325cce7SMatthew G Knepley         }
50934bf4ff8Smarkadams4         nnz1 += ncols - newcols - 1;
5109566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
5119566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
5124325cce7SMatthew G Knepley       }
5139566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5149566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5154325cce7SMatthew G Knepley     }
5169566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
5179566063dSJacob Faibussowitsch     PetscCall(PetscFree2(newCols, newVals));
5189566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */
51934bf4ff8Smarkadams4     if (nnz0 > 0) PetscCall(PetscInfo(NULL, "Filtering left %g%% edges in graph\n", 100 * (double)nnz1 / (double)nnz0));
5205afa97f7SPierre Jolivet     else PetscCall(PetscInfo(NULL, "Warning: %" PetscInt_FMT " edges to filter with %" PetscInt_FMT " rows\n", nnz0, maxRows));
521038df967SPierre Jolivet   }
5222ce66baaSPierre Jolivet   if (compress && A->ops->eliminatezeros) {
5232ce66baaSPierre Jolivet     Mat       B;
52424c7e8a4SPierre Jolivet     PetscBool flg;
5252ce66baaSPierre Jolivet 
52624c7e8a4SPierre Jolivet     PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJHIPSPARSE, MATMPIAIJHIPSPARSE, ""));
52724c7e8a4SPierre Jolivet     if (!flg) {
5282ce66baaSPierre Jolivet       PetscCall(MatEliminateZeros(A, keep));
5292ce66baaSPierre Jolivet       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
5302ce66baaSPierre Jolivet       PetscCall(MatHeaderReplace(A, &B));
5312ce66baaSPierre Jolivet     }
53224c7e8a4SPierre Jolivet   }
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5344325cce7SMatthew G Knepley }
535