16f79c3a4SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 36f79c3a4SBarry Smith 49371c9d4SSatish Balay static PetscErrorCode MatTransposeAXPY_Private(Mat Y, PetscScalar a, Mat X, MatStructure str, Mat T) { 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) { 13*11a5261eSBarry Smith PetscCall(PetscInfo(NULL, "Explicitly transposing X of type MATTRANSPOSE to perform MatAXPY()\n")); 149566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &F)); 15d54c938cSPierre Jolivet A = Y; 16d54c938cSPierre Jolivet } else { 17*11a5261eSBarry Smith PetscCall(PetscInfo(NULL, "Transposing X because Y of type MATTRANSPOSE 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) { 23*11a5261eSBarry Smith PetscCall(PetscInfo(NULL, "Explicitly Hermitian transposing X of type MATTRANSPOSE to perform MatAXPY()\n")); 249566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &F)); 25d54c938cSPierre Jolivet A = Y; 26d54c938cSPierre Jolivet } else { 27*11a5261eSBarry Smith PetscCall(PetscInfo(NULL, "Hermitian transposing X because Y of type MATTRANSPOSE 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)); 33d54c938cSPierre Jolivet PetscFunctionReturn(0); 34d54c938cSPierre Jolivet } 35d54c938cSPierre Jolivet 3606be10caSBarry Smith /*@ 3721c89e3eSBarry Smith MatAXPY - Computes Y = a*X + Y. 386f79c3a4SBarry Smith 393f9fe445SBarry Smith Logically Collective on Mat 40fee21e36SBarry Smith 4198a79cdbSBarry Smith Input Parameters: 42607cd303SBarry Smith + a - the scalar multiplier 43607cd303SBarry Smith . X - the first matrix 44607cd303SBarry Smith . Y - the second matrix 45531d0c6aSBarry 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) 4698a79cdbSBarry Smith 47*11a5261eSBarry Smith Note: 48*11a5261eSBarry Smith No operation is performed when a is zero. 49edf9a46cSStefano Zampini 502860a424SLois Curfman McInnes Level: intermediate 512860a424SLois Curfman McInnes 52db781477SPatrick Sanan .seealso: `MatAYPX()` 5306be10caSBarry Smith @*/ 549371c9d4SSatish Balay PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str) { 55646531bbSStefano Zampini PetscInt M1, M2, N1, N2; 56c1ac3661SBarry Smith PetscInt m1, m2, n1, n2; 57a683ea4aSPierre Jolivet MatType t1, t2; 58a683ea4aSPierre Jolivet PetscBool sametype, transpose; 596f79c3a4SBarry Smith 603a40ed3dSBarry Smith PetscFunctionBegin; 610700a824SBarry Smith PetscValidHeaderSpecific(Y, MAT_CLASSID, 1); 62c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y, a, 2); 63646531bbSStefano Zampini PetscValidHeaderSpecific(X, MAT_CLASSID, 3); 64646531bbSStefano Zampini PetscCheckSameComm(Y, 1, X, 3); 659566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, &M1, &N1)); 669566063dSJacob Faibussowitsch PetscCall(MatGetSize(Y, &M2, &N2)); 679566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m1, &n1)); 689566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Y, &m2, &n2)); 692c71b3e2SJacob 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); 702c71b3e2SJacob 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); 712c71b3e2SJacob Faibussowitsch PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)"); 722c71b3e2SJacob Faibussowitsch PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)"); 73edf9a46cSStefano Zampini if (a == (PetscScalar)0.0) PetscFunctionReturn(0); 74edf9a46cSStefano Zampini if (Y == X) { 759566063dSJacob Faibussowitsch PetscCall(MatScale(Y, 1.0 + a)); 76edf9a46cSStefano Zampini PetscFunctionReturn(0); 77edf9a46cSStefano Zampini } 789566063dSJacob Faibussowitsch PetscCall(MatGetType(X, &t1)); 799566063dSJacob Faibussowitsch PetscCall(MatGetType(Y, &t2)); 809566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(t1, t2, &sametype)); 819566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0)); 8293c87e65SStefano Zampini if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) { 83dbbe0bcdSBarry Smith PetscUseTypeMethod(Y, axpy, a, X, str); 84d4bb536fSBarry Smith } else { 85*11a5261eSBarry Smith PetscCall(PetscStrcmp(t1, MATTRANSPOSE, &transpose)); 86a683ea4aSPierre Jolivet if (transpose) { 879566063dSJacob Faibussowitsch PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X)); 88a683ea4aSPierre Jolivet } else { 89*11a5261eSBarry Smith PetscCall(PetscStrcmp(t2, MATTRANSPOSE, &transpose)); 90a683ea4aSPierre Jolivet if (transpose) { 919566063dSJacob Faibussowitsch PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y)); 92a683ea4aSPierre Jolivet } else { 939566063dSJacob Faibussowitsch PetscCall(MatAXPY_Basic(Y, a, X, str)); 94607cd303SBarry Smith } 95a683ea4aSPierre Jolivet } 96a683ea4aSPierre Jolivet } 979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0)); 98607cd303SBarry Smith PetscFunctionReturn(0); 99607cd303SBarry Smith } 100607cd303SBarry Smith 1019371c9d4SSatish Balay PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) { 10227afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL; 103646531bbSStefano Zampini 104646531bbSStefano Zampini PetscFunctionBegin; 105646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall)); 107646531bbSStefano Zampini if (preall) { 1089566063dSJacob Faibussowitsch PetscCall((*preall)(Y, X, B)); 109646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 110646531bbSStefano Zampini Mat preallocator; 111646531bbSStefano Zampini PetscInt r, rstart, rend; 112646531bbSStefano Zampini PetscInt m, n, M, N; 113646531bbSStefano Zampini 1149566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(Y)); 1159566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(X)); 1169566063dSJacob Faibussowitsch PetscCall(MatGetSize(Y, &M, &N)); 1179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Y, &m, &n)); 1189566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator)); 1199566063dSJacob Faibussowitsch PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 1209566063dSJacob Faibussowitsch PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap)); 1219566063dSJacob Faibussowitsch PetscCall(MatSetUp(preallocator)); 1229566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend)); 123646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 124646531bbSStefano Zampini PetscInt ncols; 125646531bbSStefano Zampini const PetscInt *row; 126646531bbSStefano Zampini const PetscScalar *vals; 127646531bbSStefano Zampini 1289566063dSJacob Faibussowitsch PetscCall(MatGetRow(Y, r, &ncols, &row, &vals)); 1299566063dSJacob Faibussowitsch PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES)); 1309566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals)); 1319566063dSJacob Faibussowitsch PetscCall(MatGetRow(X, r, &ncols, &row, &vals)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES)); 1339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals)); 134646531bbSStefano Zampini } 1359566063dSJacob Faibussowitsch PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 1369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 1379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 1389566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(Y)); 1399566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(X)); 140646531bbSStefano Zampini 1419566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B)); 1429566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name)); 1439566063dSJacob Faibussowitsch PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap)); 1449566063dSJacob Faibussowitsch PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B)); 1459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&preallocator)); 146646531bbSStefano Zampini } 147646531bbSStefano Zampini PetscFunctionReturn(0); 148646531bbSStefano Zampini } 149646531bbSStefano Zampini 1509371c9d4SSatish Balay PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str) { 151be705e3aSPierre Jolivet PetscBool isshell, isdense, isnest; 152edf9a46cSStefano Zampini 153edf9a46cSStefano Zampini PetscFunctionBegin; 1549566063dSJacob Faibussowitsch PetscCall(MatIsShell(Y, &isshell)); 155edf9a46cSStefano Zampini if (isshell) { /* MatShell has special support for AXPY */ 156edf9a46cSStefano Zampini PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure); 157edf9a46cSStefano Zampini 1589566063dSJacob Faibussowitsch PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f)); 159edf9a46cSStefano Zampini if (f) { 1609566063dSJacob Faibussowitsch PetscCall((*f)(Y, a, X, str)); 161edf9a46cSStefano Zampini PetscFunctionReturn(0); 162edf9a46cSStefano Zampini } 163edf9a46cSStefano Zampini } 16493c87e65SStefano Zampini /* no need to preallocate if Y is dense */ 1659566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, "")); 166be705e3aSPierre Jolivet if (isdense) { 1679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &isnest)); 168be705e3aSPierre Jolivet if (isnest) { 1699566063dSJacob Faibussowitsch PetscCall(MatAXPY_Dense_Nest(Y, a, X)); 170be705e3aSPierre Jolivet PetscFunctionReturn(0); 171be705e3aSPierre Jolivet } 172be705e3aSPierre Jolivet if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN; 173be705e3aSPierre Jolivet } 174a6ab3590SBarry Smith if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) { 175edf9a46cSStefano Zampini PetscInt i, start, end, j, ncols, m, n; 17638baddfdSBarry Smith const PetscInt *row; 177b3cc6726SBarry Smith PetscScalar *val; 178b3cc6726SBarry Smith const PetscScalar *vals; 179607cd303SBarry Smith 1809566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, &m, &n)); 1819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(X, &start, &end)); 1829566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(X)); 183f4df32b1SMatthew Knepley if (a == 1.0) { 184d4bb536fSBarry Smith for (i = start; i < end; i++) { 1859566063dSJacob Faibussowitsch PetscCall(MatGetRow(X, i, &ncols, &row, &vals)); 1869566063dSJacob Faibussowitsch PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES)); 1879566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals)); 188d4bb536fSBarry Smith } 189d4bb536fSBarry Smith } else { 19021a3365eSStefano Zampini PetscInt vs = 100; 19121a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 1929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(vs, &val)); 19306be10caSBarry Smith for (i = start; i < end; i++) { 1949566063dSJacob Faibussowitsch PetscCall(MatGetRow(X, i, &ncols, &row, &vals)); 19521a3365eSStefano Zampini if (vs < ncols) { 19621a3365eSStefano Zampini vs = PetscMin(2 * ncols, n); 1979566063dSJacob Faibussowitsch PetscCall(PetscRealloc(vs * sizeof(*val), &val)); 1986f79c3a4SBarry Smith } 19921a3365eSStefano Zampini for (j = 0; j < ncols; j++) val[j] = a * vals[j]; 2009566063dSJacob Faibussowitsch PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES)); 2019566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals)); 2026f79c3a4SBarry Smith } 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(val)); 204d4bb536fSBarry Smith } 2059566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(X)); 2069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY)); 2079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY)); 208edf9a46cSStefano Zampini } else { 209edf9a46cSStefano Zampini Mat B; 210edf9a46cSStefano Zampini 2119566063dSJacob Faibussowitsch PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B)); 2129566063dSJacob Faibussowitsch PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str)); 2139566063dSJacob Faibussowitsch PetscCall(MatHeaderMerge(Y, &B)); 214edf9a46cSStefano Zampini } 2153a40ed3dSBarry Smith PetscFunctionReturn(0); 2166f79c3a4SBarry Smith } 217052efed2SBarry Smith 2189371c9d4SSatish Balay PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str) { 219ec7775f6SShri Abhyankar PetscInt i, start, end, j, ncols, m, n; 220ec7775f6SShri Abhyankar const PetscInt *row; 221ec7775f6SShri Abhyankar PetscScalar *val; 222ec7775f6SShri Abhyankar const PetscScalar *vals; 223ec7775f6SShri Abhyankar 224ec7775f6SShri Abhyankar PetscFunctionBegin; 2259566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, &m, &n)); 2269566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(X, &start, &end)); 2279566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(Y)); 2289566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(X)); 229ec7775f6SShri Abhyankar if (a == 1.0) { 230ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 2319566063dSJacob Faibussowitsch PetscCall(MatGetRow(Y, i, &ncols, &row, &vals)); 2329566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES)); 2339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals)); 234ec7775f6SShri Abhyankar 2359566063dSJacob Faibussowitsch PetscCall(MatGetRow(X, i, &ncols, &row, &vals)); 2369566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES)); 2379566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals)); 238ec7775f6SShri Abhyankar } 239ec7775f6SShri Abhyankar } else { 24021a3365eSStefano Zampini PetscInt vs = 100; 24121a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 2429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(vs, &val)); 243ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 2449566063dSJacob Faibussowitsch PetscCall(MatGetRow(Y, i, &ncols, &row, &vals)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES)); 2469566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals)); 247ec7775f6SShri Abhyankar 2489566063dSJacob Faibussowitsch PetscCall(MatGetRow(X, i, &ncols, &row, &vals)); 24921a3365eSStefano Zampini if (vs < ncols) { 25021a3365eSStefano Zampini vs = PetscMin(2 * ncols, n); 2519566063dSJacob Faibussowitsch PetscCall(PetscRealloc(vs * sizeof(*val), &val)); 252ec7775f6SShri Abhyankar } 25321a3365eSStefano Zampini for (j = 0; j < ncols; j++) val[j] = a * vals[j]; 2549566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES)); 2559566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals)); 256ec7775f6SShri Abhyankar } 2579566063dSJacob Faibussowitsch PetscCall(PetscFree(val)); 258ec7775f6SShri Abhyankar } 2599566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(Y)); 2609566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(X)); 2619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 263ec7775f6SShri Abhyankar PetscFunctionReturn(0); 264ec7775f6SShri Abhyankar } 265ec7775f6SShri Abhyankar 266052efed2SBarry Smith /*@ 26787828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 268052efed2SBarry Smith 2693f9fe445SBarry Smith Neighbor-wise Collective on Mat 270fee21e36SBarry Smith 27198a79cdbSBarry Smith Input Parameters: 27298a79cdbSBarry Smith + Y - the matrices 27387828ca2SBarry Smith - a - the PetscScalar 27498a79cdbSBarry Smith 2752860a424SLois Curfman McInnes Level: intermediate 2762860a424SLois Curfman McInnes 27795452b02SPatrick Sanan Notes: 27884e19e3eSJunchao Zhang If Y is a rectangular matrix, the shift is done on the main diagonal Y_{ii} of the matrix (https://en.wikipedia.org/wiki/Main_diagonal) 27984e19e3eSJunchao Zhang 28095452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2816f33a894SBarry 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 282edf9a46cSStefano Zampini entry). No operation is performed when a is zero. 2836f33a894SBarry Smith 2840c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2850c0fd78eSBarry Smith 286db781477SPatrick Sanan .seealso: `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()` 287052efed2SBarry Smith @*/ 2889371c9d4SSatish Balay PetscErrorCode MatShift(Mat Y, PetscScalar a) { 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); 2941b71c1a7SBarry Smith if (a == 0.0) PetscFunctionReturn(0); 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)); 3003a40ed3dSBarry Smith PetscFunctionReturn(0); 301052efed2SBarry Smith } 3026d84be18SBarry Smith 3039371c9d4SSatish Balay PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is) { 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)); 31409f38230SBarry Smith PetscFunctionReturn(0); 31509f38230SBarry Smith } 31609f38230SBarry Smith 3176d84be18SBarry Smith /*@ 318f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 319f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 320f56f2b3fSBarry Smith INSERT_VALUES. 3216d84be18SBarry Smith 32248e586daSJose E. Roman Neighbor-wise Collective on Mat 32348e586daSJose E. Roman 3246d84be18SBarry Smith Input Parameters: 32598a79cdbSBarry Smith + Y - the input matrix 326f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 327f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 3286d84be18SBarry Smith 329*11a5261eSBarry Smith Note: 33095452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 3316f33a894SBarry 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 3326f33a894SBarry Smith entry). 3336f33a894SBarry Smith 3342860a424SLois Curfman McInnes Level: intermediate 3352860a424SLois Curfman McInnes 336db781477SPatrick Sanan .seealso: `MatShift()`, `MatScale()`, `MatDiagonalScale()` 3376d84be18SBarry Smith @*/ 3389371c9d4SSatish Balay PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is) { 33967576b8bSBarry Smith PetscInt matlocal, veclocal; 3406d84be18SBarry Smith 3413a40ed3dSBarry Smith PetscFunctionBegin; 3420700a824SBarry Smith PetscValidHeaderSpecific(Y, MAT_CLASSID, 1); 3430700a824SBarry Smith PetscValidHeaderSpecific(D, VEC_CLASSID, 2); 3449566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Y, &matlocal, NULL)); 3459566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(D, &veclocal)); 3465f80ce2aSJacob 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); 347dbbe0bcdSBarry Smith if (Y->ops->diagonalset) PetscUseTypeMethod(Y, diagonalset, D, is); 348dbbe0bcdSBarry Smith else PetscCall(MatDiagonalSet_Default(Y, D, is)); 3499566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 3503a40ed3dSBarry Smith PetscFunctionReturn(0); 3516d84be18SBarry Smith } 352d4bb536fSBarry Smith 353d4bb536fSBarry Smith /*@ 35404aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 355d4bb536fSBarry Smith 3563f9fe445SBarry Smith Logically on Mat 357fee21e36SBarry Smith 35898a79cdbSBarry Smith Input Parameters: 35904aac2b0SHong Zhang + a - the PetscScalar multiplier 36004aac2b0SHong Zhang . Y - the first matrix 36104aac2b0SHong Zhang . X - the second matrix 362531d0c6aSBarry 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) 363d4bb536fSBarry Smith 3642860a424SLois Curfman McInnes Level: intermediate 3652860a424SLois Curfman McInnes 366db781477SPatrick Sanan .seealso: `MatAXPY()` 367d4bb536fSBarry Smith @*/ 3689371c9d4SSatish Balay PetscErrorCode MatAYPX(Mat Y, PetscScalar a, Mat X, MatStructure str) { 3693a40ed3dSBarry Smith PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(MatScale(Y, a)); 3719566063dSJacob Faibussowitsch PetscCall(MatAXPY(Y, 1.0, X, str)); 3723a40ed3dSBarry Smith PetscFunctionReturn(0); 373d4bb536fSBarry Smith } 374b0a32e0cSBarry Smith 375b0a32e0cSBarry Smith /*@ 3760bacdadaSStefano Zampini MatComputeOperator - Computes the explicit matrix 377b0a32e0cSBarry Smith 378b0a32e0cSBarry Smith Collective on Mat 379b0a32e0cSBarry Smith 380d8d19677SJose E. Roman Input Parameters: 381186a3e20SStefano Zampini + inmat - the matrix 382186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator 383b0a32e0cSBarry Smith 384b0a32e0cSBarry Smith Output Parameter: 385a5b23f4aSJose E. Roman . mat - the explicit operator 386b0a32e0cSBarry Smith 387*11a5261eSBarry Smith Note: 388186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 389186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 390186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 391b0a32e0cSBarry Smith 392b0a32e0cSBarry Smith Level: advanced 393b0a32e0cSBarry Smith 394b0a32e0cSBarry Smith @*/ 3959371c9d4SSatish Balay PetscErrorCode MatComputeOperator(Mat inmat, MatType mattype, Mat *mat) { 396b0a32e0cSBarry Smith PetscFunctionBegin; 3970700a824SBarry Smith PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1); 398186a3e20SStefano Zampini PetscValidPointer(mat, 3); 3999566063dSJacob Faibussowitsch PetscCall(MatConvert_Shell(inmat, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat)); 40024f910e3SHong Zhang PetscFunctionReturn(0); 40124f910e3SHong Zhang } 4024325cce7SMatthew G Knepley 4034325cce7SMatthew G Knepley /*@ 4040bacdadaSStefano Zampini MatComputeOperatorTranspose - Computes the explicit matrix representation of 405f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 406f3b1f45cSBarry Smith 407f3b1f45cSBarry Smith Collective on Mat 408f3b1f45cSBarry Smith 4096b867d5aSJose E. Roman Input Parameters: 4106b867d5aSJose E. Roman + inmat - the matrix 4116b867d5aSJose E. Roman - mattype - the matrix type for the explicit operator 412f3b1f45cSBarry Smith 413f3b1f45cSBarry Smith Output Parameter: 414a5b23f4aSJose E. Roman . mat - the explicit operator transposed 415f3b1f45cSBarry Smith 416*11a5261eSBarry Smith Note: 417186a3e20SStefano Zampini This computation is done by applying the transpose of the operator to columns of the identity matrix. 418186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 419186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 420f3b1f45cSBarry Smith 421f3b1f45cSBarry Smith Level: advanced 422f3b1f45cSBarry Smith 423f3b1f45cSBarry Smith @*/ 4249371c9d4SSatish Balay PetscErrorCode MatComputeOperatorTranspose(Mat inmat, MatType mattype, Mat *mat) { 425186a3e20SStefano Zampini Mat A; 426f3b1f45cSBarry Smith 427f3b1f45cSBarry Smith PetscFunctionBegin; 428f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1); 429186a3e20SStefano Zampini PetscValidPointer(mat, 3); 4309566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(inmat, &A)); 4319566063dSJacob Faibussowitsch PetscCall(MatConvert_Shell(A, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat)); 4329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 433f3b1f45cSBarry Smith PetscFunctionReturn(0); 434f3b1f45cSBarry Smith } 435f3b1f45cSBarry Smith 436f3b1f45cSBarry Smith /*@ 4374325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4384325cce7SMatthew G Knepley 4394325cce7SMatthew G Knepley Input Parameters: 4404325cce7SMatthew G Knepley + A - The matrix 4414325cce7SMatthew G Knepley - tol - The zero tolerance 4424325cce7SMatthew G Knepley 4434325cce7SMatthew G Knepley Output Parameters: 4444325cce7SMatthew G Knepley . A - The chopped matrix 4454325cce7SMatthew G Knepley 4464325cce7SMatthew G Knepley Level: intermediate 4474325cce7SMatthew G Knepley 448db781477SPatrick Sanan .seealso: `MatCreate()`, `MatZeroEntries()` 4493fc99919SSatish Balay @*/ 4509371c9d4SSatish Balay PetscErrorCode MatChop(Mat A, PetscReal tol) { 451038df967SPierre Jolivet Mat a; 4524325cce7SMatthew G Knepley PetscScalar *newVals; 453cf5d3338SPierre Jolivet PetscInt *newCols, rStart, rEnd, numRows, maxRows, r, colMax = 0; 454038df967SPierre Jolivet PetscBool flg; 4554325cce7SMatthew G Knepley 4564325cce7SMatthew G Knepley PetscFunctionBegin; 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, "")); 458038df967SPierre Jolivet if (flg) { 4599566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(A, &a)); 4609566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(a, &r)); 4619566063dSJacob Faibussowitsch PetscCall(MatGetSize(a, &rStart, &rEnd)); 4629566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(a, &newVals)); 463038df967SPierre Jolivet for (; colMax < rEnd; ++colMax) { 464ad540459SPierre Jolivet for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) < tol ? 0.0 : newVals[maxRows + colMax * r]; 465038df967SPierre Jolivet } 4669566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(a, &newVals)); 467038df967SPierre Jolivet } else { 4689566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 4699566063dSJacob Faibussowitsch PetscCall(MatGetRowUpperTriangular(A)); 4704325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 4714325cce7SMatthew G Knepley PetscInt ncols; 4724325cce7SMatthew G Knepley 4739566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, NULL, NULL)); 4745f80ce2aSJacob Faibussowitsch colMax = PetscMax(colMax, ncols); 4759566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL)); 4764325cce7SMatthew G Knepley } 4774325cce7SMatthew G Knepley numRows = rEnd - rStart; 4781c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A))); 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(colMax, &newCols, colMax, &newVals)); 4809566063dSJacob Faibussowitsch PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */ 4819566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 482cf5d3338SPierre Jolivet /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd() */ 483cf5d3338SPierre Jolivet /* that are potentially called many times depending on the distribution of A */ 4844325cce7SMatthew G Knepley for (r = rStart; r < rStart + maxRows; ++r) { 4854325cce7SMatthew G Knepley const PetscScalar *vals; 4864325cce7SMatthew G Knepley const PetscInt *cols; 487fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 4884325cce7SMatthew G Knepley 4894325cce7SMatthew G Knepley if (r < rEnd) { 4909566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, &cols, &vals)); 4914325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 4924325cce7SMatthew G Knepley newCols[c] = cols[c]; 4934325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 4944325cce7SMatthew G Knepley } 495fad45679SMatthew G. Knepley newcols = ncols; 4969566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals)); 4979566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES)); 4984325cce7SMatthew G Knepley } 4999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 5014325cce7SMatthew G Knepley } 5029566063dSJacob Faibussowitsch PetscCall(MatRestoreRowUpperTriangular(A)); 5039566063dSJacob Faibussowitsch PetscCall(PetscFree2(newCols, newVals)); 5049566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */ 505038df967SPierre Jolivet } 5064325cce7SMatthew G Knepley PetscFunctionReturn(0); 5074325cce7SMatthew G Knepley } 508