xref: /petsc/src/snes/tutorials/ex6.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4421ceaSFande Kong 
2c4421ceaSFande Kong static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3c4421ceaSFande Kong This example employs a user-defined reasonview routine.\n\n";
4c4421ceaSFande Kong 
5c4421ceaSFande Kong #include <petscsnes.h>
6c4421ceaSFande Kong 
7c4421ceaSFande Kong /*
8c4421ceaSFande Kong    User-defined routines
9c4421ceaSFande Kong */
10c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
11c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
12c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormInitialGuess(Vec);
13c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MySNESConvergedReasonView(SNES, void *);
14c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MyKSPConvergedReasonView(KSP, void *);
15c4421ceaSFande Kong 
16c4421ceaSFande Kong /*
17c4421ceaSFande Kong    User-defined context for monitoring
18c4421ceaSFande Kong */
19c4421ceaSFande Kong typedef struct {
20c4421ceaSFande Kong   PetscViewer viewer;
21c4421ceaSFande Kong } ReasonViewCtx;
22c4421ceaSFande Kong 
23d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
24d71ae5a4SJacob Faibussowitsch {
25c4421ceaSFande Kong   SNES          snes;       /* SNES context */
26c4421ceaSFande Kong   KSP           ksp;        /* KSP context */
27c4421ceaSFande Kong   Vec           x, r, F, U; /* vectors */
28c4421ceaSFande Kong   Mat           J;          /* Jacobian matrix */
29c4421ceaSFande Kong   ReasonViewCtx monP;       /* monitoring context */
30c4421ceaSFande Kong   PetscInt      its, n = 5, i;
31c4421ceaSFande Kong   PetscMPIInt   size;
32c4421ceaSFande Kong   PetscScalar   h, xp, v;
33c4421ceaSFande Kong   MPI_Comm      comm;
34c4421ceaSFande Kong 
35327415f7SBarry Smith   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
38be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
40c4421ceaSFande Kong   h    = 1.0 / (n - 1);
41c4421ceaSFande Kong   comm = PETSC_COMM_WORLD;
42c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43c4421ceaSFande Kong      Create nonlinear solver context
44c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45c4421ceaSFande Kong 
469566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
47c4421ceaSFande Kong 
48c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49c4421ceaSFande Kong      Create vector data structures; set function evaluation routine
50c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51c4421ceaSFande Kong   /*
52c4421ceaSFande Kong      Note that we form 1 vector from scratch and then duplicate as needed.
53c4421ceaSFande Kong   */
549566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &x));
559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, n));
569566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &F));
599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
60c4421ceaSFande Kong 
61c4421ceaSFande Kong   /*
62c4421ceaSFande Kong      Set function evaluation routine and vector
63c4421ceaSFande Kong   */
649566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F));
65c4421ceaSFande Kong 
66c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4421ceaSFande Kong      Create matrix data structure; set Jacobian evaluation routine
68c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69c4421ceaSFande Kong 
709566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &J));
719566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, n, n));
729566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
739566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
74c4421ceaSFande Kong 
75c4421ceaSFande Kong   /*
76c4421ceaSFande Kong      Set Jacobian matrix data structure and default Jacobian evaluation
77c4421ceaSFande Kong      routine. User can override with:
78c4421ceaSFande Kong      -snes_fd : default finite differencing approximation of Jacobian
79c4421ceaSFande Kong      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
80c4421ceaSFande Kong                 (unless user explicitly sets preconditioner)
81c4421ceaSFande Kong      -snes_mf_operator : form preconditioning matrix as set by the user,
82c4421ceaSFande Kong                          but use matrix-free approx for Jacobian-vector
83c4421ceaSFande Kong                          products within Newton-Krylov method
84c4421ceaSFande Kong   */
85c4421ceaSFande Kong 
869566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
87c4421ceaSFande Kong 
88c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4421ceaSFande Kong      Customize nonlinear solver; set runtime options
90c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91c4421ceaSFande Kong 
92c4421ceaSFande Kong   /*
93c4421ceaSFande Kong      Set an optional user-defined reasonview routine
94c4421ceaSFande Kong   */
959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(comm, &monP.viewer));
96d5b43468SJose E. Roman   /* Just make sure we can not repeat adding the same function
97c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
98c4421ceaSFande Kong    */
9948a46eb9SPierre Jolivet   for (i = 0; i < 4; i++) PetscCall(SNESConvergedReasonViewSet(snes, MySNESConvergedReasonView, &monP, 0));
1009566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
101d5b43468SJose E. Roman   /* Just make sure we can not repeat adding the same function
102c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
103c4421ceaSFande Kong    */
10448a46eb9SPierre Jolivet   for (i = 0; i < 4; i++) PetscCall(KSPConvergedReasonViewSet(ksp, MyKSPConvergedReasonView, &monP, 0));
105c4421ceaSFande Kong   /*
106c4421ceaSFande Kong      Set SNES/KSP/KSP/PC runtime options, e.g.,
107c4421ceaSFande Kong          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
108c4421ceaSFande Kong   */
1099566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
110c4421ceaSFande Kong 
111c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112c4421ceaSFande Kong      Initialize application:
113c4421ceaSFande Kong      Store right-hand-side of PDE and exact solution
114c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115c4421ceaSFande Kong 
116c4421ceaSFande Kong   xp = 0.0;
117c4421ceaSFande Kong   for (i = 0; i < n; i++) {
118c4421ceaSFande Kong     v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
1199566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES));
120c4421ceaSFande Kong     v = xp * xp * xp;
1219566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES));
122c4421ceaSFande Kong     xp += h;
123c4421ceaSFande Kong   }
124c4421ceaSFande Kong 
125c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4421ceaSFande Kong      Evaluate initial guess; then solve nonlinear system
127c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4421ceaSFande Kong   /*
129c4421ceaSFande Kong      Note: The user should initialize the vector, x, with the initial guess
130c4421ceaSFande Kong      for the nonlinear solver prior to calling SNESSolve().  In particular,
131c4421ceaSFande Kong      to employ an initial guess of zero, the user should explicitly set
132c4421ceaSFande Kong      this vector to zero by calling VecSet().
133c4421ceaSFande Kong   */
1349566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1359566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1369566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
137c4421ceaSFande Kong 
138c4421ceaSFande Kong   /*
139c4421ceaSFande Kong      Free work space.  All PETSc objects should be destroyed when they
140c4421ceaSFande Kong      are no longer needed.
141c4421ceaSFande Kong   */
1429371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1439371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1449371c9d4SSatish Balay   PetscCall(VecDestroy(&U));
1459371c9d4SSatish Balay   PetscCall(VecDestroy(&F));
1469371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1479371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
149b122ec5aSJacob Faibussowitsch   return 0;
150c4421ceaSFande Kong }
151f6dfbefdSBarry Smith 
152c4421ceaSFande Kong /*
153c4421ceaSFande Kong    FormInitialGuess - Computes initial guess.
154c4421ceaSFande Kong 
155c4421ceaSFande Kong    Input/Output Parameter:
156c4421ceaSFande Kong .  x - the solution vector
157c4421ceaSFande Kong */
158d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
159d71ae5a4SJacob Faibussowitsch {
160c4421ceaSFande Kong   PetscScalar pfive = .50;
161*3ba16761SJacob Faibussowitsch 
162*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1639566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
164*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165c4421ceaSFande Kong }
166f6dfbefdSBarry Smith 
167c4421ceaSFande Kong /*
168c4421ceaSFande Kong    FormFunction - Evaluates nonlinear function, F(x).
169c4421ceaSFande Kong 
170c4421ceaSFande Kong    Input Parameters:
171c4421ceaSFande Kong .  snes - the SNES context
172c4421ceaSFande Kong .  x - input vector
173c4421ceaSFande Kong .  ctx - optional user-defined context, as set by SNESSetFunction()
174c4421ceaSFande Kong 
175c4421ceaSFande Kong    Output Parameter:
176c4421ceaSFande Kong .  f - function vector
177c4421ceaSFande Kong 
178c4421ceaSFande Kong    Note:
179c4421ceaSFande Kong    The user-defined context can contain any application-specific data
180c4421ceaSFande Kong    needed for the function evaluation (such as various parameters, work
181c4421ceaSFande Kong    vectors, and grid information).  In this program the context is just
182c4421ceaSFande Kong    a vector containing the right-hand-side of the discretized PDE.
183c4421ceaSFande Kong  */
184c4421ceaSFande Kong 
185d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx)
186d71ae5a4SJacob Faibussowitsch {
187c4421ceaSFande Kong   Vec                g = (Vec)ctx;
188c4421ceaSFande Kong   const PetscScalar *xx, *gg;
189c4421ceaSFande Kong   PetscScalar       *ff, d;
190c4421ceaSFande Kong   PetscInt           i, n;
191c4421ceaSFande Kong 
192*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
193c4421ceaSFande Kong   /*
194c4421ceaSFande Kong      Get pointers to vector data.
195c4421ceaSFande Kong        - For default PETSc vectors, VecGetArray() returns a pointer to
196c4421ceaSFande Kong          the data array.  Otherwise, the routine is implementation dependent.
197c4421ceaSFande Kong        - You MUST call VecRestoreArray() when you no longer need access to
198c4421ceaSFande Kong          the array.
199c4421ceaSFande Kong   */
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
2029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g, &gg));
203c4421ceaSFande Kong 
204c4421ceaSFande Kong   /*
205c4421ceaSFande Kong      Compute function
206c4421ceaSFande Kong   */
2079566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2089371c9d4SSatish Balay   d     = (PetscReal)(n - 1);
2099371c9d4SSatish Balay   d     = d * d;
210c4421ceaSFande Kong   ff[0] = xx[0];
211c4421ceaSFande Kong   for (i = 1; i < n - 1; i++) ff[i] = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - gg[i];
212c4421ceaSFande Kong   ff[n - 1] = xx[n - 1] - 1.0;
213c4421ceaSFande Kong 
214c4421ceaSFande Kong   /*
215c4421ceaSFande Kong      Restore vectors
216c4421ceaSFande Kong   */
2179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
2199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g, &gg));
220*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221c4421ceaSFande Kong }
222c4421ceaSFande Kong /* ------------------------------------------------------------------- */
223c4421ceaSFande Kong /*
224c4421ceaSFande Kong    FormJacobian - Evaluates Jacobian matrix.
225c4421ceaSFande Kong 
226c4421ceaSFande Kong    Input Parameters:
227c4421ceaSFande Kong .  snes - the SNES context
228c4421ceaSFande Kong .  x - input vector
229c4421ceaSFande Kong .  dummy - optional user-defined context (not used here)
230c4421ceaSFande Kong 
231c4421ceaSFande Kong    Output Parameters:
232c4421ceaSFande Kong .  jac - Jacobian matrix
233c4421ceaSFande Kong .  B - optionally different preconditioning matrix
234c4421ceaSFande Kong 
235c4421ceaSFande Kong */
236c4421ceaSFande Kong 
237d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
238d71ae5a4SJacob Faibussowitsch {
239c4421ceaSFande Kong   const PetscScalar *xx;
240c4421ceaSFande Kong   PetscScalar        A[3], d;
241c4421ceaSFande Kong   PetscInt           i, n, j[3];
242c4421ceaSFande Kong 
243*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
244c4421ceaSFande Kong   /*
245c4421ceaSFande Kong      Get pointer to vector data
246c4421ceaSFande Kong   */
2479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
248c4421ceaSFande Kong 
249c4421ceaSFande Kong   /*
250c4421ceaSFande Kong      Compute Jacobian entries and insert into matrix.
251c4421ceaSFande Kong       - Note that in this case we set all elements for a particular
252c4421ceaSFande Kong         row at once.
253c4421ceaSFande Kong   */
2549566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2559371c9d4SSatish Balay   d = (PetscReal)(n - 1);
2569371c9d4SSatish Balay   d = d * d;
257c4421ceaSFande Kong 
258c4421ceaSFande Kong   /*
259c4421ceaSFande Kong      Interior grid points
260c4421ceaSFande Kong   */
261c4421ceaSFande Kong   for (i = 1; i < n - 1; i++) {
2629371c9d4SSatish Balay     j[0] = i - 1;
2639371c9d4SSatish Balay     j[1] = i;
2649371c9d4SSatish Balay     j[2] = i + 1;
2659371c9d4SSatish Balay     A[0] = A[2] = d;
2669371c9d4SSatish Balay     A[1]        = -2.0 * d + 2.0 * xx[i];
2679566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
268c4421ceaSFande Kong   }
269c4421ceaSFande Kong 
270c4421ceaSFande Kong   /*
271c4421ceaSFande Kong      Boundary points
272c4421ceaSFande Kong   */
2739371c9d4SSatish Balay   i    = 0;
2749371c9d4SSatish Balay   A[0] = 1.0;
275c4421ceaSFande Kong 
2769566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
277c4421ceaSFande Kong 
2789371c9d4SSatish Balay   i    = n - 1;
2799371c9d4SSatish Balay   A[0] = 1.0;
280c4421ceaSFande Kong 
2819566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
282c4421ceaSFande Kong 
283c4421ceaSFande Kong   /*
284c4421ceaSFande Kong      Restore vector
285c4421ceaSFande Kong   */
2869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
287c4421ceaSFande Kong 
288c4421ceaSFande Kong   /*
289c4421ceaSFande Kong      Assemble matrix
290c4421ceaSFande Kong   */
2919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
293c4421ceaSFande Kong   if (jac != B) {
2949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2959566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
296c4421ceaSFande Kong   }
297*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
298c4421ceaSFande Kong }
299c4421ceaSFande Kong 
300d71ae5a4SJacob Faibussowitsch PetscErrorCode MySNESConvergedReasonView(SNES snes, void *ctx)
301d71ae5a4SJacob Faibussowitsch {
302c4421ceaSFande Kong   ReasonViewCtx      *monP   = (ReasonViewCtx *)ctx;
303c4421ceaSFande Kong   PetscViewer         viewer = monP->viewer;
304c4421ceaSFande Kong   SNESConvergedReason reason;
305c4421ceaSFande Kong   const char         *strreason;
306c4421ceaSFande Kong 
307*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
3089566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes, &reason));
3099566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReasonString(snes, &strreason));
3109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Customized SNES converged reason view\n"));
3119566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, 1));
312c4421ceaSFande Kong   if (reason > 0) {
3139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Converged due to %s\n", strreason));
314c4421ceaSFande Kong   } else if (reason <= 0) {
3159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Did not converge due to %s\n", strreason));
316c4421ceaSFande Kong   }
3179566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, 1));
318*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
319c4421ceaSFande Kong }
320c4421ceaSFande Kong 
321d71ae5a4SJacob Faibussowitsch PetscErrorCode MyKSPConvergedReasonView(KSP ksp, void *ctx)
322d71ae5a4SJacob Faibussowitsch {
323c4421ceaSFande Kong   ReasonViewCtx     *monP   = (ReasonViewCtx *)ctx;
324c4421ceaSFande Kong   PetscViewer        viewer = monP->viewer;
325c4421ceaSFande Kong   KSPConvergedReason reason;
326c4421ceaSFande Kong   const char        *reasonstr;
327c4421ceaSFande Kong 
328*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
3299566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(ksp, &reason));
3309566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReasonString(ksp, &reasonstr));
3319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, 2));
3329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Customized KSP converged reason view\n"));
3339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, 1));
334c4421ceaSFande Kong   if (reason > 0) {
3359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Converged due to %s\n", reasonstr));
336c4421ceaSFande Kong   } else if (reason <= 0) {
3379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Did not converge due to %s\n", reasonstr));
338c4421ceaSFande Kong   }
3399566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, 3));
340*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
341c4421ceaSFande Kong }
342c4421ceaSFande Kong 
343c4421ceaSFande Kong /*TEST
344c4421ceaSFande Kong 
345c4421ceaSFande Kong    test:
346c4421ceaSFande Kong       suffix: 1
347c4421ceaSFande Kong       nsize: 1
348bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
349c4421ceaSFande Kong 
350c4421ceaSFande Kong    test:
351c4421ceaSFande Kong       suffix: 2
352c4421ceaSFande Kong       nsize: 1
353c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel
354bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
355c4421ceaSFande Kong 
356c4421ceaSFande Kong    test:
357c4421ceaSFande Kong       suffix: 3
358c4421ceaSFande Kong       nsize: 1
359c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel -ksp_converged_reason
360bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
361c4421ceaSFande Kong 
362c4421ceaSFande Kong    test:
363c4421ceaSFande Kong       suffix: 4
364c4421ceaSFande Kong       nsize: 1
365c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel
366bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
367c4421ceaSFande Kong 
368c4421ceaSFande Kong    test:
369c4421ceaSFande Kong       suffix: 5
370c4421ceaSFande Kong       nsize: 1
371c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel -snes_converged_reason
372bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
373c4421ceaSFande Kong 
374c4421ceaSFande Kong TEST*/
375