xref: /petsc/src/mat/utils/axpy.c (revision 2e37cc59afff36ccb5a1b1fe6479c31382db32d5)
1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
26f79c3a4SBarry Smith 
3d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTransposeAXPY_Private(Mat Y, PetscScalar a, Mat X, MatStructure str, Mat T)
4d71ae5a4SJacob Faibussowitsch {
5d54c938cSPierre Jolivet   Mat         A, F;
64166f22eSPierre Jolivet   PetscScalar vshift, vscale;
75f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat, Mat *);
8d54c938cSPierre Jolivet 
9d54c938cSPierre Jolivet   PetscFunctionBegin;
104166f22eSPierre Jolivet   if (T == X) PetscCall(MatShellGetScalingShifts(T, &vshift, &vscale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
114166f22eSPierre Jolivet   else {
124166f22eSPierre Jolivet     vshift = 0.0;
134166f22eSPierre Jolivet     vscale = 1.0;
144166f22eSPierre Jolivet   }
159566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)T, "MatTransposeGetMat_C", &f));
16d54c938cSPierre Jolivet   if (f) {
179566063dSJacob Faibussowitsch     PetscCall(MatTransposeGetMat(T, &A));
18d54c938cSPierre Jolivet     if (T == X) {
19013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Explicitly transposing X of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
209566063dSJacob Faibussowitsch       PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &F));
21d54c938cSPierre Jolivet       A = Y;
22d54c938cSPierre Jolivet     } else {
23013e2dc7SBarry Smith       PetscCall(PetscInfo(NULL, "Transposing X because Y of type MATTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
249566063dSJacob Faibussowitsch       PetscCall(MatTranspose(X, MAT_INITIAL_MATRIX, &F));
25d54c938cSPierre Jolivet     }
26d54c938cSPierre Jolivet   } else {
279566063dSJacob Faibussowitsch     PetscCall(MatHermitianTransposeGetMat(T, &A));
28d54c938cSPierre Jolivet     if (T == X) {
2952cd20c4SPierre Jolivet       PetscCall(PetscInfo(NULL, "Explicitly Hermitian transposing X of type MATHERMITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
309566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &F));
31d54c938cSPierre Jolivet       A = Y;
32d54c938cSPierre Jolivet     } else {
3352cd20c4SPierre Jolivet       PetscCall(PetscInfo(NULL, "Hermitian transposing X because Y of type MATHERMITIANTRANSPOSEVIRTUAL to perform MatAXPY()\n"));
349566063dSJacob Faibussowitsch       PetscCall(MatHermitianTranspose(X, MAT_INITIAL_MATRIX, &F));
35d54c938cSPierre Jolivet     }
36d54c938cSPierre Jolivet   }
374166f22eSPierre Jolivet   PetscCall(MatAXPY(A, a * vscale, F, str));
384166f22eSPierre Jolivet   PetscCall(MatShift(A, a * vshift));
399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41d54c938cSPierre Jolivet }
42d54c938cSPierre Jolivet 
4306be10caSBarry Smith /*@
4421c89e3eSBarry Smith   MatAXPY - Computes Y = a*X + Y.
456f79c3a4SBarry Smith 
46c3339decSBarry Smith   Logically Collective
47fee21e36SBarry Smith 
4898a79cdbSBarry Smith   Input Parameters:
49607cd303SBarry Smith + a   - the scalar multiplier
50607cd303SBarry Smith . X   - the first matrix
51607cd303SBarry Smith . Y   - the second matrix
5227430b45SBarry Smith - str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN`, `UNKNOWN_NONZERO_PATTERN`, or `SUBSET_NONZERO_PATTERN` (nonzeros of `X` is a subset of `Y`'s)
53edf9a46cSStefano Zampini 
542860a424SLois Curfman McInnes   Level: intermediate
552860a424SLois Curfman McInnes 
561cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAYPX()`
5706be10caSBarry Smith  @*/
58d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str)
59d71ae5a4SJacob Faibussowitsch {
60646531bbSStefano Zampini   PetscInt  M1, M2, N1, N2;
61c1ac3661SBarry Smith   PetscInt  m1, m2, n1, n2;
62*2e37cc59SJose E. Roman   PetscBool sametype, transpose, nogetrow;
636f79c3a4SBarry Smith 
643a40ed3dSBarry Smith   PetscFunctionBegin;
650700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
66c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y, a, 2);
67646531bbSStefano Zampini   PetscValidHeaderSpecific(X, MAT_CLASSID, 3);
68646531bbSStefano Zampini   PetscCheckSameComm(Y, 1, X, 3);
699566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M1, &N1));
709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M2, &N2));
719566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m1, &n1));
729566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m2, &n2));
732c71b3e2SJacob Faibussowitsch   PetscCheck(M1 == M2 && N1 == N2, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_SIZ, "Non conforming matrix add: global sizes X %" PetscInt_FMT " x %" PetscInt_FMT ", Y %" PetscInt_FMT " x %" PetscInt_FMT, M1, N1, M2, N2);
742c71b3e2SJacob Faibussowitsch   PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Non conforming matrix add: local sizes X %" PetscInt_FMT " x %" PetscInt_FMT ", Y %" PetscInt_FMT " x %" PetscInt_FMT, m1, n1, m2, n2);
752c71b3e2SJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)");
762c71b3e2SJacob Faibussowitsch   PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)");
773ba16761SJacob Faibussowitsch   if (a == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
78edf9a46cSStefano Zampini   if (Y == X) {
799566063dSJacob Faibussowitsch     PetscCall(MatScale(Y, 1.0 + a));
803ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
81edf9a46cSStefano Zampini   }
82013e2dc7SBarry Smith   PetscCall(PetscObjectObjectTypeCompare((PetscObject)X, (PetscObject)Y, &sametype));
839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0));
8493c87e65SStefano Zampini   if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) {
85dbbe0bcdSBarry Smith     PetscUseTypeMethod(Y, axpy, a, X, str);
86d4bb536fSBarry Smith   } else {
87013e2dc7SBarry Smith     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
88a683ea4aSPierre Jolivet     if (transpose) {
899566063dSJacob Faibussowitsch       PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X));
90a683ea4aSPierre Jolivet     } else {
91013e2dc7SBarry Smith       PetscCall(PetscObjectTypeCompareAny((PetscObject)Y, &transpose, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
92a683ea4aSPierre Jolivet       if (transpose) {
939566063dSJacob Faibussowitsch         PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y));
94a683ea4aSPierre Jolivet       } else {
95*2e37cc59SJose E. Roman         PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &nogetrow, MATSCALAPACK, MATELEMENTAL, ""));
96*2e37cc59SJose E. Roman         if (nogetrow) { /* Avoid MatAXPY_Basic() due to missing MatGetRow() */
97*2e37cc59SJose E. Roman           Mat     C;
98*2e37cc59SJose E. Roman           MatType type;
99*2e37cc59SJose E. Roman 
100*2e37cc59SJose E. Roman           PetscCall(MatGetType(Y, &type));
101*2e37cc59SJose E. Roman           PetscCall(MatConvert(X, type, MAT_INITIAL_MATRIX, &C));
102*2e37cc59SJose E. Roman           PetscCall(MatAXPY(Y, a, C, str));
103*2e37cc59SJose E. Roman           PetscCall(MatDestroy(&C));
104*2e37cc59SJose E. Roman         } else {
1059566063dSJacob Faibussowitsch           PetscCall(MatAXPY_Basic(Y, a, X, str));
106607cd303SBarry Smith         }
107a683ea4aSPierre Jolivet       }
108a683ea4aSPierre Jolivet     }
109*2e37cc59SJose E. Roman   }
1109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0));
1113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112607cd303SBarry Smith }
113607cd303SBarry Smith 
114d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
115d71ae5a4SJacob Faibussowitsch {
11627afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL;
117646531bbSStefano Zampini 
118646531bbSStefano Zampini   PetscFunctionBegin;
119646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
1209566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall));
121646531bbSStefano Zampini   if (preall) {
1229566063dSJacob Faibussowitsch     PetscCall((*preall)(Y, X, B));
123646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
124646531bbSStefano Zampini     Mat      preallocator;
125646531bbSStefano Zampini     PetscInt r, rstart, rend;
126646531bbSStefano Zampini     PetscInt m, n, M, N;
127646531bbSStefano Zampini 
1289566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
1299566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
1309566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Y, &M, &N));
1319566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Y, &m, &n));
1329566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator));
1339566063dSJacob Faibussowitsch     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1349566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap));
1359566063dSJacob Faibussowitsch     PetscCall(MatSetUp(preallocator));
1369566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
137646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
138646531bbSStefano Zampini       PetscInt           ncols;
139646531bbSStefano Zampini       const PetscInt    *row;
140646531bbSStefano Zampini       const PetscScalar *vals;
141646531bbSStefano Zampini 
1429566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, r, &ncols, &row, &vals));
1439566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1449566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals));
1459566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, r, &ncols, &row, &vals));
1469566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
1479566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals));
148646531bbSStefano Zampini     }
1499566063dSJacob Faibussowitsch     PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1509566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1519566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
1529566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
1539566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
154646531bbSStefano Zampini 
1559566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B));
1569566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name));
1579566063dSJacob Faibussowitsch     PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap));
1589566063dSJacob Faibussowitsch     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B));
1599566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&preallocator));
160646531bbSStefano Zampini   }
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162646531bbSStefano Zampini }
163646531bbSStefano Zampini 
164d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str)
165d71ae5a4SJacob Faibussowitsch {
166be705e3aSPierre Jolivet   PetscBool isshell, isdense, isnest;
167edf9a46cSStefano Zampini 
168edf9a46cSStefano Zampini   PetscFunctionBegin;
1699566063dSJacob Faibussowitsch   PetscCall(MatIsShell(Y, &isshell));
170edf9a46cSStefano Zampini   if (isshell) { /* MatShell has special support for AXPY */
171edf9a46cSStefano Zampini     PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure);
172edf9a46cSStefano Zampini 
1739566063dSJacob Faibussowitsch     PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f));
174edf9a46cSStefano Zampini     if (f) {
1759566063dSJacob Faibussowitsch       PetscCall((*f)(Y, a, X, str));
1763ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
177edf9a46cSStefano Zampini     }
178edf9a46cSStefano Zampini   }
17993c87e65SStefano Zampini   /* no need to preallocate if Y is dense */
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, ""));
181be705e3aSPierre Jolivet   if (isdense) {
1829566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &isnest));
183be705e3aSPierre Jolivet     if (isnest) {
1849566063dSJacob Faibussowitsch       PetscCall(MatAXPY_Dense_Nest(Y, a, X));
1853ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
186be705e3aSPierre Jolivet     }
187be705e3aSPierre Jolivet     if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN;
188be705e3aSPierre Jolivet   }
189a6ab3590SBarry Smith   if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) {
190edf9a46cSStefano Zampini     PetscInt           i, start, end, j, ncols, m, n;
19138baddfdSBarry Smith     const PetscInt    *row;
192b3cc6726SBarry Smith     PetscScalar       *val;
193b3cc6726SBarry Smith     const PetscScalar *vals;
194607cd303SBarry Smith 
1959566063dSJacob Faibussowitsch     PetscCall(MatGetSize(X, &m, &n));
1969566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(X, &start, &end));
1979566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
198f4df32b1SMatthew Knepley     if (a == 1.0) {
199d4bb536fSBarry Smith       for (i = start; i < end; i++) {
2009566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2019566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES));
2029566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
203d4bb536fSBarry Smith       }
204d4bb536fSBarry Smith     } else {
20521a3365eSStefano Zampini       PetscInt vs = 100;
20621a3365eSStefano Zampini       /* realloc if needed, as this function may be used in parallel */
2079566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(vs, &val));
20806be10caSBarry Smith       for (i = start; i < end; i++) {
2099566063dSJacob Faibussowitsch         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
21021a3365eSStefano Zampini         if (vs < ncols) {
21121a3365eSStefano Zampini           vs = PetscMin(2 * ncols, n);
2129566063dSJacob Faibussowitsch           PetscCall(PetscRealloc(vs * sizeof(*val), &val));
2136f79c3a4SBarry Smith         }
21421a3365eSStefano Zampini         for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2159566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES));
2169566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
2176f79c3a4SBarry Smith       }
2189566063dSJacob Faibussowitsch       PetscCall(PetscFree(val));
219d4bb536fSBarry Smith     }
2209566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
2219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
2229566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
223edf9a46cSStefano Zampini   } else {
224edf9a46cSStefano Zampini     Mat B;
225edf9a46cSStefano Zampini 
2269566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B));
2279566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
2289566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y, &B));
229edf9a46cSStefano Zampini   }
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2316f79c3a4SBarry Smith }
232052efed2SBarry Smith 
233d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str)
234d71ae5a4SJacob Faibussowitsch {
235ec7775f6SShri Abhyankar   PetscInt           i, start, end, j, ncols, m, n;
236ec7775f6SShri Abhyankar   const PetscInt    *row;
237ec7775f6SShri Abhyankar   PetscScalar       *val;
238ec7775f6SShri Abhyankar   const PetscScalar *vals;
239ec7775f6SShri Abhyankar 
240ec7775f6SShri Abhyankar   PetscFunctionBegin;
2419566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &m, &n));
2429566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(X, &start, &end));
2439566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(Y));
2449566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(X));
245ec7775f6SShri Abhyankar   if (a == 1.0) {
246ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2479566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2489566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2499566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
250ec7775f6SShri Abhyankar 
2519566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
2529566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2539566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
254ec7775f6SShri Abhyankar     }
255ec7775f6SShri Abhyankar   } else {
25621a3365eSStefano Zampini     PetscInt vs = 100;
25721a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
2589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(vs, &val));
259ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
2609566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
2619566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
2629566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
263ec7775f6SShri Abhyankar 
2649566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
26521a3365eSStefano Zampini       if (vs < ncols) {
26621a3365eSStefano Zampini         vs = PetscMin(2 * ncols, n);
2679566063dSJacob Faibussowitsch         PetscCall(PetscRealloc(vs * sizeof(*val), &val));
268ec7775f6SShri Abhyankar       }
26921a3365eSStefano Zampini       for (j = 0; j < ncols; j++) val[j] = a * vals[j];
2709566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES));
2719566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
272ec7775f6SShri Abhyankar     }
2739566063dSJacob Faibussowitsch     PetscCall(PetscFree(val));
274ec7775f6SShri Abhyankar   }
2759566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(Y));
2769566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(X));
2779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
280ec7775f6SShri Abhyankar }
281ec7775f6SShri Abhyankar 
282052efed2SBarry Smith /*@
2832ef1f0ffSBarry Smith   MatShift - Computes `Y =  Y + a I`, where `a` is a `PetscScalar`
284052efed2SBarry Smith 
285c3339decSBarry Smith   Neighbor-wise Collective
286fee21e36SBarry Smith 
28798a79cdbSBarry Smith   Input Parameters:
28827430b45SBarry Smith + Y - the matrix
28927430b45SBarry Smith - a - the `PetscScalar`
29098a79cdbSBarry Smith 
2912860a424SLois Curfman McInnes   Level: intermediate
2922860a424SLois Curfman McInnes 
29395452b02SPatrick Sanan   Notes:
29427430b45SBarry 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)
29584e19e3eSJunchao Zhang 
2962ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2976f33a894SBarry 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
298edf9a46cSStefano Zampini   entry). No operation is performed when a is zero.
2996f33a894SBarry Smith 
30027430b45SBarry Smith   To form Y = Y + diag(V) use `MatDiagonalSet()`
3010c0fd78eSBarry Smith 
3021cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()`
303052efed2SBarry Smith  @*/
304d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift(Mat Y, PetscScalar a)
305d71ae5a4SJacob Faibussowitsch {
3063a40ed3dSBarry Smith   PetscFunctionBegin;
3070700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3085f80ce2aSJacob Faibussowitsch   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3095f80ce2aSJacob Faibussowitsch   PetscCheck(!Y->factortype, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3104994cf47SJed Brown   MatCheckPreallocated(Y, 1);
3113ba16761SJacob Faibussowitsch   if (a == 0.0) PetscFunctionReturn(PETSC_SUCCESS);
312b50b34bdSBarry Smith 
313dbbe0bcdSBarry Smith   if (Y->ops->shift) PetscUseTypeMethod(Y, shift, a);
3149566063dSJacob Faibussowitsch   else PetscCall(MatShift_Basic(Y, a));
3157d68702bSBarry Smith 
3169566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318052efed2SBarry Smith }
3196d84be18SBarry Smith 
320d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is)
321d71ae5a4SJacob Faibussowitsch {
32267576b8bSBarry Smith   PetscInt           i, start, end;
323fa2e578bSStefano Zampini   const PetscScalar *v;
32409f38230SBarry Smith 
32509f38230SBarry Smith   PetscFunctionBegin;
3269566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y, &start, &end));
3279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(D, &v));
32848a46eb9SPierre Jolivet   for (i = start; i < end; i++) PetscCall(MatSetValues(Y, 1, &i, 1, &i, v + i - start, is));
3299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(D, &v));
3309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
3319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33309f38230SBarry Smith }
33409f38230SBarry Smith 
3356d84be18SBarry Smith /*@
3362ef1f0ffSBarry Smith   MatDiagonalSet - Computes `Y` = `Y` + `D`, where `D` is a diagonal matrix
3372ef1f0ffSBarry Smith   that is represented as a vector. Or Y[i,i] = D[i] if `InsertMode` is
3382ef1f0ffSBarry Smith   `INSERT_VALUES`.
3396d84be18SBarry Smith 
340c3339decSBarry Smith   Neighbor-wise Collective
34148e586daSJose E. Roman 
3426d84be18SBarry Smith   Input Parameters:
34398a79cdbSBarry Smith + Y  - the input matrix
344f56f2b3fSBarry Smith . D  - the diagonal matrix, represented as a vector
345fe59aa6dSJacob Faibussowitsch - is - `INSERT_VALUES` or `ADD_VALUES`
3466f33a894SBarry Smith 
3472860a424SLois Curfman McInnes   Level: intermediate
3482860a424SLois Curfman McInnes 
3492ef1f0ffSBarry Smith   Note:
3502ef1f0ffSBarry Smith   If the matrix `Y` is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3512ef1f0ffSBarry 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
3522ef1f0ffSBarry Smith   entry).
3532ef1f0ffSBarry Smith 
3541cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatShift()`, `MatScale()`, `MatDiagonalScale()`
3556d84be18SBarry Smith @*/
356d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is)
357d71ae5a4SJacob Faibussowitsch {
35867576b8bSBarry Smith   PetscInt matlocal, veclocal;
3596d84be18SBarry Smith 
3603a40ed3dSBarry Smith   PetscFunctionBegin;
3610700a824SBarry Smith   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
3620700a824SBarry Smith   PetscValidHeaderSpecific(D, VEC_CLASSID, 2);
3632fbb122aSJose E. Roman   MatCheckPreallocated(Y, 1);
3649566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &matlocal, NULL));
3659566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(D, &veclocal));
3665f80ce2aSJacob 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);
367dbbe0bcdSBarry Smith   if (Y->ops->diagonalset) PetscUseTypeMethod(Y, diagonalset, D, is);
368dbbe0bcdSBarry Smith   else PetscCall(MatDiagonalSet_Default(Y, D, is));
3699566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3716d84be18SBarry Smith }
372d4bb536fSBarry Smith 
373d4bb536fSBarry Smith /*@
37404aac2b0SHong Zhang   MatAYPX - Computes Y = a*Y + X.
375d4bb536fSBarry Smith 
3762ef1f0ffSBarry Smith   Logically Collective
377fee21e36SBarry Smith 
37898a79cdbSBarry Smith   Input Parameters:
3792ef1f0ffSBarry Smith + a   - the `PetscScalar` multiplier
38004aac2b0SHong Zhang . Y   - the first matrix
38104aac2b0SHong Zhang . X   - the second matrix
382aaa8cc7dSPierre 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)
383d4bb536fSBarry Smith 
3842860a424SLois Curfman McInnes   Level: intermediate
3852860a424SLois Curfman McInnes 
3861cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAXPY()`
387d4bb536fSBarry Smith  @*/
388d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAYPX(Mat Y, PetscScalar a, Mat X, MatStructure str)
389d71ae5a4SJacob Faibussowitsch {
3903a40ed3dSBarry Smith   PetscFunctionBegin;
3919566063dSJacob Faibussowitsch   PetscCall(MatScale(Y, a));
3929566063dSJacob Faibussowitsch   PetscCall(MatAXPY(Y, 1.0, X, str));
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
394d4bb536fSBarry Smith }
395b0a32e0cSBarry Smith 
3965d83a8b1SBarry Smith /*@
3970bacdadaSStefano Zampini   MatComputeOperator - Computes the explicit matrix
398b0a32e0cSBarry Smith 
399c3339decSBarry Smith   Collective
400b0a32e0cSBarry Smith 
401d8d19677SJose E. Roman   Input Parameters:
402186a3e20SStefano Zampini + inmat   - the matrix
403186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator
404b0a32e0cSBarry Smith 
405b0a32e0cSBarry Smith   Output Parameter:
406a5b23f4aSJose E. Roman . mat - the explicit  operator
407b0a32e0cSBarry Smith 
4082ef1f0ffSBarry Smith   Level: advanced
4092ef1f0ffSBarry Smith 
41011a5261eSBarry Smith   Note:
411f13dfd9eSBarry Smith   This computation is done by applying the operator to columns of the identity matrix.
412186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
4132ef1f0ffSBarry Smith   Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
414b0a32e0cSBarry Smith 
4151cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperatorTranspose()`
416b0a32e0cSBarry Smith @*/
417d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperator(Mat inmat, MatType mattype, Mat *mat)
418d71ae5a4SJacob Faibussowitsch {
419b0a32e0cSBarry Smith   PetscFunctionBegin;
4200700a824SBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
4214f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
4229566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(inmat, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42424f910e3SHong Zhang }
4254325cce7SMatthew G Knepley 
4265d83a8b1SBarry Smith /*@
4270bacdadaSStefano Zampini   MatComputeOperatorTranspose - Computes the explicit matrix representation of
4282ef1f0ffSBarry Smith   a give matrix that can apply `MatMultTranspose()`
429f3b1f45cSBarry Smith 
430c3339decSBarry Smith   Collective
431f3b1f45cSBarry Smith 
4326b867d5aSJose E. Roman   Input Parameters:
4336b867d5aSJose E. Roman + inmat   - the matrix
4346b867d5aSJose E. Roman - mattype - the matrix type for the explicit operator
435f3b1f45cSBarry Smith 
436f3b1f45cSBarry Smith   Output Parameter:
437a5b23f4aSJose E. Roman . mat - the explicit  operator transposed
438f3b1f45cSBarry Smith 
4392ef1f0ffSBarry Smith   Level: advanced
4402ef1f0ffSBarry Smith 
44111a5261eSBarry Smith   Note:
442186a3e20SStefano Zampini   This computation is done by applying the transpose of the operator to columns of the identity matrix.
443186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
4442ef1f0ffSBarry Smith   Currently, this routine uses a dense matrix format if `mattype` == `NULL`.
445f3b1f45cSBarry Smith 
4461cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatMult()`, `MatComputeOperator()`
447f3b1f45cSBarry Smith @*/
448d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeOperatorTranspose(Mat inmat, MatType mattype, Mat *mat)
449d71ae5a4SJacob Faibussowitsch {
450186a3e20SStefano Zampini   Mat A;
451f3b1f45cSBarry Smith 
452f3b1f45cSBarry Smith   PetscFunctionBegin;
453f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
4544f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
4559566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(inmat, &A));
4569566063dSJacob Faibussowitsch   PetscCall(MatConvert_Shell(A, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
4579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
459f3b1f45cSBarry Smith }
460f3b1f45cSBarry Smith 
461f3b1f45cSBarry Smith /*@
4622ce66baaSPierre 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
4634325cce7SMatthew G Knepley 
4644325cce7SMatthew G Knepley   Input Parameters:
4654325cce7SMatthew G Knepley + A        - The matrix
4662ce66baaSPierre Jolivet . tol      - The zero tolerance
4672ce66baaSPierre 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
4682ce66baaSPierre 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
4694325cce7SMatthew G Knepley 
4704325cce7SMatthew G Knepley   Level: intermediate
4714325cce7SMatthew G Knepley 
472d9b70fcdSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatZeroEntries()`, `MatEliminateZeros()`, `VecFilter()`
4733fc99919SSatish Balay  @*/
4742ce66baaSPierre Jolivet PetscErrorCode MatFilter(Mat A, PetscReal tol, PetscBool compress, PetscBool keep)
475d71ae5a4SJacob Faibussowitsch {
476038df967SPierre Jolivet   Mat          a;
4774325cce7SMatthew G Knepley   PetscScalar *newVals;
47834bf4ff8Smarkadams4   PetscInt    *newCols, rStart, rEnd, maxRows, r, colMax = 0, nnz0 = 0, nnz1 = 0;
479038df967SPierre Jolivet   PetscBool    flg;
4804325cce7SMatthew G Knepley 
4814325cce7SMatthew G Knepley   PetscFunctionBegin;
4829566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, ""));
483038df967SPierre Jolivet   if (flg) {
4849566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLocalMatrix(A, &a));
4859566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(a, &r));
4869566063dSJacob Faibussowitsch     PetscCall(MatGetSize(a, &rStart, &rEnd));
4879566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(a, &newVals));
488038df967SPierre Jolivet     for (; colMax < rEnd; ++colMax) {
48982fa6941SPierre Jolivet       for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) <= tol ? 0.0 : newVals[maxRows + colMax * r];
490038df967SPierre Jolivet     }
4919566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(a, &newVals));
492038df967SPierre Jolivet   } else {
4933da7c217SPierre Jolivet     const PetscInt *ranges;
4943da7c217SPierre Jolivet     PetscMPIInt     rank, size;
4953da7c217SPierre Jolivet 
4963da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
4973da7c217SPierre Jolivet     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
4983da7c217SPierre Jolivet     PetscCall(MatGetOwnershipRanges(A, &ranges));
4993da7c217SPierre Jolivet     rStart = ranges[rank];
5003da7c217SPierre Jolivet     rEnd   = ranges[rank + 1];
5019566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
5024325cce7SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
5034325cce7SMatthew G Knepley       PetscInt ncols;
5044325cce7SMatthew G Knepley 
5059566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, r, &ncols, NULL, NULL));
5065f80ce2aSJacob Faibussowitsch       colMax = PetscMax(colMax, ncols);
5079566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL));
5084325cce7SMatthew G Knepley     }
5093da7c217SPierre Jolivet     maxRows = 0;
5103da7c217SPierre Jolivet     for (r = 0; r < size; ++r) maxRows = PetscMax(maxRows, ranges[r + 1] - ranges[r]);
5113da7c217SPierre Jolivet     PetscCall(PetscCalloc2(colMax, &newCols, colMax, &newVals));
5129566063dSJacob Faibussowitsch     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */
5139566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
514cf5d3338SPierre Jolivet     /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd()             */
515cf5d3338SPierre Jolivet     /* that are potentially called many times depending on the distribution of A */
5164325cce7SMatthew G Knepley     for (r = rStart; r < rStart + maxRows; ++r) {
5173da7c217SPierre Jolivet       if (r < rEnd) {
5184325cce7SMatthew G Knepley         const PetscScalar *vals;
5194325cce7SMatthew G Knepley         const PetscInt    *cols;
5203da7c217SPierre Jolivet         PetscInt           ncols, newcols = 0, c;
5214325cce7SMatthew G Knepley 
5229566063dSJacob Faibussowitsch         PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
52334bf4ff8Smarkadams4         nnz0 += ncols - 1;
5244325cce7SMatthew G Knepley         for (c = 0; c < ncols; ++c) {
52582fa6941SPierre Jolivet           if (PetscUnlikely(PetscAbsScalar(vals[c]) <= tol)) newCols[newcols++] = cols[c];
5264325cce7SMatthew G Knepley         }
52734bf4ff8Smarkadams4         nnz1 += ncols - newcols - 1;
5289566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
5299566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
5304325cce7SMatthew G Knepley       }
5319566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5329566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5334325cce7SMatthew G Knepley     }
5349566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
5359566063dSJacob Faibussowitsch     PetscCall(PetscFree2(newCols, newVals));
5369566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */
53734bf4ff8Smarkadams4     if (nnz0 > 0) PetscCall(PetscInfo(NULL, "Filtering left %g%% edges in graph\n", 100 * (double)nnz1 / (double)nnz0));
5385afa97f7SPierre Jolivet     else PetscCall(PetscInfo(NULL, "Warning: %" PetscInt_FMT " edges to filter with %" PetscInt_FMT " rows\n", nnz0, maxRows));
539038df967SPierre Jolivet   }
5402ce66baaSPierre Jolivet   if (compress && A->ops->eliminatezeros) {
5412ce66baaSPierre Jolivet     Mat       B;
54224c7e8a4SPierre Jolivet     PetscBool flg;
5432ce66baaSPierre Jolivet 
54424c7e8a4SPierre Jolivet     PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJHIPSPARSE, MATMPIAIJHIPSPARSE, ""));
54524c7e8a4SPierre Jolivet     if (!flg) {
5462ce66baaSPierre Jolivet       PetscCall(MatEliminateZeros(A, keep));
5472ce66baaSPierre Jolivet       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
5482ce66baaSPierre Jolivet       PetscCall(MatHeaderReplace(A, &B));
5492ce66baaSPierre Jolivet     }
55024c7e8a4SPierre Jolivet   }
5513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5524325cce7SMatthew G Knepley }
553