xref: /petsc/src/mat/utils/axpy.c (revision 11a5261e40035b7c793f2783a2ba6c7cd4f3b077)
1 
2 #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
3 
4 static PetscErrorCode MatTransposeAXPY_Private(Mat Y, PetscScalar a, Mat X, MatStructure str, Mat T) {
5   Mat A, F;
6   PetscErrorCode (*f)(Mat, Mat *);
7 
8   PetscFunctionBegin;
9   PetscCall(PetscObjectQueryFunction((PetscObject)T, "MatTransposeGetMat_C", &f));
10   if (f) {
11     PetscCall(MatTransposeGetMat(T, &A));
12     if (T == X) {
13       PetscCall(PetscInfo(NULL, "Explicitly transposing X of type MATTRANSPOSE to perform MatAXPY()\n"));
14       PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &F));
15       A = Y;
16     } else {
17       PetscCall(PetscInfo(NULL, "Transposing X because Y of type MATTRANSPOSE to perform MatAXPY()\n"));
18       PetscCall(MatTranspose(X, MAT_INITIAL_MATRIX, &F));
19     }
20   } else {
21     PetscCall(MatHermitianTransposeGetMat(T, &A));
22     if (T == X) {
23       PetscCall(PetscInfo(NULL, "Explicitly Hermitian transposing X of type MATTRANSPOSE to perform MatAXPY()\n"));
24       PetscCall(MatHermitianTranspose(A, MAT_INITIAL_MATRIX, &F));
25       A = Y;
26     } else {
27       PetscCall(PetscInfo(NULL, "Hermitian transposing X because Y of type MATTRANSPOSE to perform MatAXPY()\n"));
28       PetscCall(MatHermitianTranspose(X, MAT_INITIAL_MATRIX, &F));
29     }
30   }
31   PetscCall(MatAXPY(A, a, F, str));
32   PetscCall(MatDestroy(&F));
33   PetscFunctionReturn(0);
34 }
35 
36 /*@
37    MatAXPY - Computes Y = a*X + Y.
38 
39    Logically  Collective on Mat
40 
41    Input Parameters:
42 +  a - the scalar multiplier
43 .  X - the first matrix
44 .  Y - the second matrix
45 -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, UNKNOWN_NONZERO_PATTERN, or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)
46 
47    Note:
48    No operation is performed when a is zero.
49 
50    Level: intermediate
51 
52 .seealso: `MatAYPX()`
53  @*/
54 PetscErrorCode MatAXPY(Mat Y, PetscScalar a, Mat X, MatStructure str) {
55   PetscInt  M1, M2, N1, N2;
56   PetscInt  m1, m2, n1, n2;
57   MatType   t1, t2;
58   PetscBool sametype, transpose;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
62   PetscValidLogicalCollectiveScalar(Y, a, 2);
63   PetscValidHeaderSpecific(X, MAT_CLASSID, 3);
64   PetscCheckSameComm(Y, 1, X, 3);
65   PetscCall(MatGetSize(X, &M1, &N1));
66   PetscCall(MatGetSize(Y, &M2, &N2));
67   PetscCall(MatGetLocalSize(X, &m1, &n1));
68   PetscCall(MatGetLocalSize(Y, &m2, &n2));
69   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);
70   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);
71   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (Y)");
72   PetscCheck(X->assembled, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix (X)");
73   if (a == (PetscScalar)0.0) PetscFunctionReturn(0);
74   if (Y == X) {
75     PetscCall(MatScale(Y, 1.0 + a));
76     PetscFunctionReturn(0);
77   }
78   PetscCall(MatGetType(X, &t1));
79   PetscCall(MatGetType(Y, &t2));
80   PetscCall(PetscStrcmp(t1, t2, &sametype));
81   PetscCall(PetscLogEventBegin(MAT_AXPY, Y, 0, 0, 0));
82   if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) {
83     PetscUseTypeMethod(Y, axpy, a, X, str);
84   } else {
85     PetscCall(PetscStrcmp(t1, MATTRANSPOSE, &transpose));
86     if (transpose) {
87       PetscCall(MatTransposeAXPY_Private(Y, a, X, str, X));
88     } else {
89       PetscCall(PetscStrcmp(t2, MATTRANSPOSE, &transpose));
90       if (transpose) {
91         PetscCall(MatTransposeAXPY_Private(Y, a, X, str, Y));
92       } else {
93         PetscCall(MatAXPY_Basic(Y, a, X, str));
94       }
95     }
96   }
97   PetscCall(PetscLogEventEnd(MAT_AXPY, Y, 0, 0, 0));
98   PetscFunctionReturn(0);
99 }
100 
101 PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) {
102   PetscErrorCode (*preall)(Mat, Mat, Mat *) = NULL;
103 
104   PetscFunctionBegin;
105   /* look for any available faster alternative to the general preallocator */
106   PetscCall(PetscObjectQueryFunction((PetscObject)Y, "MatAXPYGetPreallocation_C", &preall));
107   if (preall) {
108     PetscCall((*preall)(Y, X, B));
109   } else { /* Use MatPrellocator, assumes same row-col distribution */
110     Mat      preallocator;
111     PetscInt r, rstart, rend;
112     PetscInt m, n, M, N;
113 
114     PetscCall(MatGetRowUpperTriangular(Y));
115     PetscCall(MatGetRowUpperTriangular(X));
116     PetscCall(MatGetSize(Y, &M, &N));
117     PetscCall(MatGetLocalSize(Y, &m, &n));
118     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &preallocator));
119     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
120     PetscCall(MatSetLayouts(preallocator, Y->rmap, Y->cmap));
121     PetscCall(MatSetUp(preallocator));
122     PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
123     for (r = rstart; r < rend; ++r) {
124       PetscInt           ncols;
125       const PetscInt    *row;
126       const PetscScalar *vals;
127 
128       PetscCall(MatGetRow(Y, r, &ncols, &row, &vals));
129       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
130       PetscCall(MatRestoreRow(Y, r, &ncols, &row, &vals));
131       PetscCall(MatGetRow(X, r, &ncols, &row, &vals));
132       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
133       PetscCall(MatRestoreRow(X, r, &ncols, &row, &vals));
134     }
135     PetscCall(MatSetOption(preallocator, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
136     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
137     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
138     PetscCall(MatRestoreRowUpperTriangular(Y));
139     PetscCall(MatRestoreRowUpperTriangular(X));
140 
141     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), B));
142     PetscCall(MatSetType(*B, ((PetscObject)Y)->type_name));
143     PetscCall(MatSetLayouts(*B, Y->rmap, Y->cmap));
144     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *B));
145     PetscCall(MatDestroy(&preallocator));
146   }
147   PetscFunctionReturn(0);
148 }
149 
150 PetscErrorCode MatAXPY_Basic(Mat Y, PetscScalar a, Mat X, MatStructure str) {
151   PetscBool isshell, isdense, isnest;
152 
153   PetscFunctionBegin;
154   PetscCall(MatIsShell(Y, &isshell));
155   if (isshell) { /* MatShell has special support for AXPY */
156     PetscErrorCode (*f)(Mat, PetscScalar, Mat, MatStructure);
157 
158     PetscCall(MatGetOperation(Y, MATOP_AXPY, (void (**)(void)) & f));
159     if (f) {
160       PetscCall((*f)(Y, a, X, str));
161       PetscFunctionReturn(0);
162     }
163   }
164   /* no need to preallocate if Y is dense */
165   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isdense, MATSEQDENSE, MATMPIDENSE, ""));
166   if (isdense) {
167     PetscCall(PetscObjectTypeCompare((PetscObject)X, MATNEST, &isnest));
168     if (isnest) {
169       PetscCall(MatAXPY_Dense_Nest(Y, a, X));
170       PetscFunctionReturn(0);
171     }
172     if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN;
173   }
174   if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) {
175     PetscInt           i, start, end, j, ncols, m, n;
176     const PetscInt    *row;
177     PetscScalar       *val;
178     const PetscScalar *vals;
179 
180     PetscCall(MatGetSize(X, &m, &n));
181     PetscCall(MatGetOwnershipRange(X, &start, &end));
182     PetscCall(MatGetRowUpperTriangular(X));
183     if (a == 1.0) {
184       for (i = start; i < end; i++) {
185         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
186         PetscCall(MatSetValues(Y, 1, &i, ncols, row, vals, ADD_VALUES));
187         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
188       }
189     } else {
190       PetscInt vs = 100;
191       /* realloc if needed, as this function may be used in parallel */
192       PetscCall(PetscMalloc1(vs, &val));
193       for (i = start; i < end; i++) {
194         PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
195         if (vs < ncols) {
196           vs = PetscMin(2 * ncols, n);
197           PetscCall(PetscRealloc(vs * sizeof(*val), &val));
198         }
199         for (j = 0; j < ncols; j++) val[j] = a * vals[j];
200         PetscCall(MatSetValues(Y, 1, &i, ncols, row, val, ADD_VALUES));
201         PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
202       }
203       PetscCall(PetscFree(val));
204     }
205     PetscCall(MatRestoreRowUpperTriangular(X));
206     PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
207     PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
208   } else {
209     Mat B;
210 
211     PetscCall(MatAXPY_Basic_Preallocate(Y, X, &B));
212     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
213     PetscCall(MatHeaderMerge(Y, &B));
214   }
215   PetscFunctionReturn(0);
216 }
217 
218 PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B, Mat Y, PetscScalar a, Mat X, MatStructure str) {
219   PetscInt           i, start, end, j, ncols, m, n;
220   const PetscInt    *row;
221   PetscScalar       *val;
222   const PetscScalar *vals;
223 
224   PetscFunctionBegin;
225   PetscCall(MatGetSize(X, &m, &n));
226   PetscCall(MatGetOwnershipRange(X, &start, &end));
227   PetscCall(MatGetRowUpperTriangular(Y));
228   PetscCall(MatGetRowUpperTriangular(X));
229   if (a == 1.0) {
230     for (i = start; i < end; i++) {
231       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
232       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
233       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
234 
235       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
236       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
237       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
238     }
239   } else {
240     PetscInt vs = 100;
241     /* realloc if needed, as this function may be used in parallel */
242     PetscCall(PetscMalloc1(vs, &val));
243     for (i = start; i < end; i++) {
244       PetscCall(MatGetRow(Y, i, &ncols, &row, &vals));
245       PetscCall(MatSetValues(B, 1, &i, ncols, row, vals, ADD_VALUES));
246       PetscCall(MatRestoreRow(Y, i, &ncols, &row, &vals));
247 
248       PetscCall(MatGetRow(X, i, &ncols, &row, &vals));
249       if (vs < ncols) {
250         vs = PetscMin(2 * ncols, n);
251         PetscCall(PetscRealloc(vs * sizeof(*val), &val));
252       }
253       for (j = 0; j < ncols; j++) val[j] = a * vals[j];
254       PetscCall(MatSetValues(B, 1, &i, ncols, row, val, ADD_VALUES));
255       PetscCall(MatRestoreRow(X, i, &ncols, &row, &vals));
256     }
257     PetscCall(PetscFree(val));
258   }
259   PetscCall(MatRestoreRowUpperTriangular(Y));
260   PetscCall(MatRestoreRowUpperTriangular(X));
261   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
262   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
263   PetscFunctionReturn(0);
264 }
265 
266 /*@
267    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
268 
269    Neighbor-wise Collective on Mat
270 
271    Input Parameters:
272 +  Y - the matrices
273 -  a - the PetscScalar
274 
275    Level: intermediate
276 
277    Notes:
278     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)
279 
280     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
281    fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an
282    entry). No operation is performed when a is zero.
283 
284    To form Y = Y + diag(V) use MatDiagonalSet()
285 
286 .seealso: `MatDiagonalSet()`, `MatScale()`, `MatDiagonalScale()`
287  @*/
288 PetscErrorCode MatShift(Mat Y, PetscScalar a) {
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
291   PetscCheck(Y->assembled, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
292   PetscCheck(!Y->factortype, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
293   MatCheckPreallocated(Y, 1);
294   if (a == 0.0) PetscFunctionReturn(0);
295 
296   if (Y->ops->shift) PetscUseTypeMethod(Y, shift, a);
297   else PetscCall(MatShift_Basic(Y, a));
298 
299   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
300   PetscFunctionReturn(0);
301 }
302 
303 PetscErrorCode MatDiagonalSet_Default(Mat Y, Vec D, InsertMode is) {
304   PetscInt           i, start, end;
305   const PetscScalar *v;
306 
307   PetscFunctionBegin;
308   PetscCall(MatGetOwnershipRange(Y, &start, &end));
309   PetscCall(VecGetArrayRead(D, &v));
310   for (i = start; i < end; i++) PetscCall(MatSetValues(Y, 1, &i, 1, &i, v + i - start, is));
311   PetscCall(VecRestoreArrayRead(D, &v));
312   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
313   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
314   PetscFunctionReturn(0);
315 }
316 
317 /*@
318    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
319    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
320    INSERT_VALUES.
321 
322    Neighbor-wise Collective on Mat
323 
324    Input Parameters:
325 +  Y - the input matrix
326 .  D - the diagonal matrix, represented as a vector
327 -  i - INSERT_VALUES or ADD_VALUES
328 
329    Note:
330     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
331    fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an
332    entry).
333 
334    Level: intermediate
335 
336 .seealso: `MatShift()`, `MatScale()`, `MatDiagonalScale()`
337 @*/
338 PetscErrorCode MatDiagonalSet(Mat Y, Vec D, InsertMode is) {
339   PetscInt matlocal, veclocal;
340 
341   PetscFunctionBegin;
342   PetscValidHeaderSpecific(Y, MAT_CLASSID, 1);
343   PetscValidHeaderSpecific(D, VEC_CLASSID, 2);
344   PetscCall(MatGetLocalSize(Y, &matlocal, NULL));
345   PetscCall(VecGetLocalSize(D, &veclocal));
346   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);
347   if (Y->ops->diagonalset) PetscUseTypeMethod(Y, diagonalset, D, is);
348   else PetscCall(MatDiagonalSet_Default(Y, D, is));
349   PetscCall(PetscObjectStateIncrease((PetscObject)Y));
350   PetscFunctionReturn(0);
351 }
352 
353 /*@
354    MatAYPX - Computes Y = a*Y + X.
355 
356    Logically on Mat
357 
358    Input Parameters:
359 +  a - the PetscScalar multiplier
360 .  Y - the first matrix
361 .  X - the second matrix
362 -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, UNKNOWN_NONZERO_PATTERN, or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)
363 
364    Level: intermediate
365 
366 .seealso: `MatAXPY()`
367  @*/
368 PetscErrorCode MatAYPX(Mat Y, PetscScalar a, Mat X, MatStructure str) {
369   PetscFunctionBegin;
370   PetscCall(MatScale(Y, a));
371   PetscCall(MatAXPY(Y, 1.0, X, str));
372   PetscFunctionReturn(0);
373 }
374 
375 /*@
376     MatComputeOperator - Computes the explicit matrix
377 
378     Collective on Mat
379 
380     Input Parameters:
381 +   inmat - the matrix
382 -   mattype - the matrix type for the explicit operator
383 
384     Output Parameter:
385 .   mat - the explicit  operator
386 
387     Note:
388     This computation is done by applying the operators to columns of the identity matrix.
389     This routine is costly in general, and is recommended for use only with relatively small systems.
390     Currently, this routine uses a dense matrix format if mattype == NULL.
391 
392     Level: advanced
393 
394 @*/
395 PetscErrorCode MatComputeOperator(Mat inmat, MatType mattype, Mat *mat) {
396   PetscFunctionBegin;
397   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
398   PetscValidPointer(mat, 3);
399   PetscCall(MatConvert_Shell(inmat, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
400   PetscFunctionReturn(0);
401 }
402 
403 /*@
404     MatComputeOperatorTranspose - Computes the explicit matrix representation of
405         a give matrix that can apply MatMultTranspose()
406 
407     Collective on Mat
408 
409     Input Parameters:
410 +   inmat - the matrix
411 -   mattype - the matrix type for the explicit operator
412 
413     Output Parameter:
414 .   mat - the explicit  operator transposed
415 
416     Note:
417     This computation is done by applying the transpose of the operator to columns of the identity matrix.
418     This routine is costly in general, and is recommended for use only with relatively small systems.
419     Currently, this routine uses a dense matrix format if mattype == NULL.
420 
421     Level: advanced
422 
423 @*/
424 PetscErrorCode MatComputeOperatorTranspose(Mat inmat, MatType mattype, Mat *mat) {
425   Mat A;
426 
427   PetscFunctionBegin;
428   PetscValidHeaderSpecific(inmat, MAT_CLASSID, 1);
429   PetscValidPointer(mat, 3);
430   PetscCall(MatCreateTranspose(inmat, &A));
431   PetscCall(MatConvert_Shell(A, mattype ? mattype : MATDENSE, MAT_INITIAL_MATRIX, mat));
432   PetscCall(MatDestroy(&A));
433   PetscFunctionReturn(0);
434 }
435 
436 /*@
437   MatChop - Set all values in the matrix less than the tolerance to zero
438 
439   Input Parameters:
440 + A   - The matrix
441 - tol - The zero tolerance
442 
443   Output Parameters:
444 . A - The chopped matrix
445 
446   Level: intermediate
447 
448 .seealso: `MatCreate()`, `MatZeroEntries()`
449  @*/
450 PetscErrorCode MatChop(Mat A, PetscReal tol) {
451   Mat          a;
452   PetscScalar *newVals;
453   PetscInt    *newCols, rStart, rEnd, numRows, maxRows, r, colMax = 0;
454   PetscBool    flg;
455 
456   PetscFunctionBegin;
457   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, ""));
458   if (flg) {
459     PetscCall(MatDenseGetLocalMatrix(A, &a));
460     PetscCall(MatDenseGetLDA(a, &r));
461     PetscCall(MatGetSize(a, &rStart, &rEnd));
462     PetscCall(MatDenseGetArray(a, &newVals));
463     for (; colMax < rEnd; ++colMax) {
464       for (maxRows = 0; maxRows < rStart; ++maxRows) newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) < tol ? 0.0 : newVals[maxRows + colMax * r];
465     }
466     PetscCall(MatDenseRestoreArray(a, &newVals));
467   } else {
468     PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
469     PetscCall(MatGetRowUpperTriangular(A));
470     for (r = rStart; r < rEnd; ++r) {
471       PetscInt ncols;
472 
473       PetscCall(MatGetRow(A, r, &ncols, NULL, NULL));
474       colMax = PetscMax(colMax, ncols);
475       PetscCall(MatRestoreRow(A, r, &ncols, NULL, NULL));
476     }
477     numRows = rEnd - rStart;
478     PetscCall(MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A)));
479     PetscCall(PetscMalloc2(colMax, &newCols, colMax, &newVals));
480     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg)); /* cache user-defined value */
481     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
482     /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd()             */
483     /* that are potentially called many times depending on the distribution of A */
484     for (r = rStart; r < rStart + maxRows; ++r) {
485       const PetscScalar *vals;
486       const PetscInt    *cols;
487       PetscInt           ncols, newcols, c;
488 
489       if (r < rEnd) {
490         PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
491         for (c = 0; c < ncols; ++c) {
492           newCols[c] = cols[c];
493           newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
494         }
495         newcols = ncols;
496         PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
497         PetscCall(MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
498       }
499       PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
500       PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
501     }
502     PetscCall(MatRestoreRowUpperTriangular(A));
503     PetscCall(PetscFree2(newCols, newVals));
504     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg)); /* reset option to its user-defined value */
505   }
506   PetscFunctionReturn(0);
507 }
508