xref: /petsc/src/snes/tutorials/ex6.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4421ceaSFande Kong static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
2c4421ceaSFande Kong This example employs a user-defined reasonview routine.\n\n";
3c4421ceaSFande Kong 
4c4421ceaSFande Kong #include <petscsnes.h>
5c4421ceaSFande Kong 
6c4421ceaSFande Kong /*
7c4421ceaSFande Kong    User-defined routines
8c4421ceaSFande Kong */
9c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
10c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
11c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormInitialGuess(Vec);
12c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MySNESConvergedReasonView(SNES, void *);
13c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MyKSPConvergedReasonView(KSP, void *);
14c4421ceaSFande Kong 
15c4421ceaSFande Kong /*
16c4421ceaSFande Kong    User-defined context for monitoring
17c4421ceaSFande Kong */
18c4421ceaSFande Kong typedef struct {
19c4421ceaSFande Kong   PetscViewer viewer;
20c4421ceaSFande Kong } ReasonViewCtx;
21c4421ceaSFande Kong 
main(int argc,char ** argv)22d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
23d71ae5a4SJacob Faibussowitsch {
24c4421ceaSFande Kong   SNES          snes;       /* SNES context */
25c4421ceaSFande Kong   KSP           ksp;        /* KSP context */
26c4421ceaSFande Kong   Vec           x, r, F, U; /* vectors */
27c4421ceaSFande Kong   Mat           J;          /* Jacobian matrix */
28c4421ceaSFande Kong   ReasonViewCtx monP;       /* monitoring context */
29c4421ceaSFande Kong   PetscInt      its, n = 5, i;
30c4421ceaSFande Kong   PetscMPIInt   size;
31c4421ceaSFande Kong   PetscScalar   h, xp, v;
32c4421ceaSFande Kong   MPI_Comm      comm;
33c4421ceaSFande Kong 
34327415f7SBarry Smith   PetscFunctionBeginUser;
35c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
37be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
39c4421ceaSFande Kong   h    = 1.0 / (n - 1);
40c4421ceaSFande Kong   comm = PETSC_COMM_WORLD;
41c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42c4421ceaSFande Kong      Create nonlinear solver context
43c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44c4421ceaSFande Kong 
459566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
46c4421ceaSFande Kong 
47c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48c4421ceaSFande Kong      Create vector data structures; set function evaluation routine
49c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50c4421ceaSFande Kong   /*
51c4421ceaSFande Kong      Note that we form 1 vector from scratch and then duplicate as needed.
52c4421ceaSFande Kong   */
539566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &x));
549566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, n));
559566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &F));
589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
59c4421ceaSFande Kong 
60c4421ceaSFande Kong   /*
61c4421ceaSFande Kong      Set function evaluation routine and vector
62c4421ceaSFande Kong   */
639566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F));
64c4421ceaSFande Kong 
65c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4421ceaSFande Kong      Create matrix data structure; set Jacobian evaluation routine
67c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68c4421ceaSFande Kong 
699566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &J));
709566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, n, n));
719566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
729566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
73c4421ceaSFande Kong 
74c4421ceaSFande Kong   /*
75c4421ceaSFande Kong      Set Jacobian matrix data structure and default Jacobian evaluation
76c4421ceaSFande Kong      routine. User can override with:
77c4421ceaSFande Kong      -snes_fd : default finite differencing approximation of Jacobian
78c4421ceaSFande Kong      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
79c4421ceaSFande Kong                 (unless user explicitly sets preconditioner)
807addb90fSBarry Smith      -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
81c4421ceaSFande Kong                          but use matrix-free approx for Jacobian-vector
82c4421ceaSFande Kong                          products within Newton-Krylov method
83c4421ceaSFande Kong   */
84c4421ceaSFande Kong 
859566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
86c4421ceaSFande Kong 
87c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88c4421ceaSFande Kong      Customize nonlinear solver; set runtime options
89c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90c4421ceaSFande Kong 
91c4421ceaSFande Kong   /*
92c4421ceaSFande Kong      Set an optional user-defined reasonview routine
93c4421ceaSFande Kong   */
949566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(comm, &monP.viewer));
95d5b43468SJose E. Roman   /* Just make sure we can not repeat adding the same function
96da81f932SPierre Jolivet    * PETSc will be able to ignore the repeated function
97c4421ceaSFande Kong    */
9848a46eb9SPierre Jolivet   for (i = 0; i < 4; i++) PetscCall(SNESConvergedReasonViewSet(snes, MySNESConvergedReasonView, &monP, 0));
999566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
100d5b43468SJose E. Roman   /* Just make sure we can not repeat adding the same function
101da81f932SPierre Jolivet    * PETSc will be able to ignore the repeated function
102c4421ceaSFande Kong    */
10348a46eb9SPierre Jolivet   for (i = 0; i < 4; i++) PetscCall(KSPConvergedReasonViewSet(ksp, MyKSPConvergedReasonView, &monP, 0));
104c4421ceaSFande Kong   /*
105c4421ceaSFande Kong      Set SNES/KSP/KSP/PC runtime options, e.g.,
106c4421ceaSFande Kong          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
107c4421ceaSFande Kong   */
1089566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
109c4421ceaSFande Kong 
110c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4421ceaSFande Kong      Initialize application:
112dd8e379bSPierre Jolivet      Store right-hand side of PDE and exact solution
113c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114c4421ceaSFande Kong 
115c4421ceaSFande Kong   xp = 0.0;
116c4421ceaSFande Kong   for (i = 0; i < n; i++) {
117c4421ceaSFande Kong     v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
1189566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES));
119c4421ceaSFande Kong     v = xp * xp * xp;
1209566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES));
121c4421ceaSFande Kong     xp += h;
122c4421ceaSFande Kong   }
123c4421ceaSFande Kong 
124c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4421ceaSFande Kong      Evaluate initial guess; then solve nonlinear system
126c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127c4421ceaSFande Kong   /*
128c4421ceaSFande Kong      Note: The user should initialize the vector, x, with the initial guess
129c4421ceaSFande Kong      for the nonlinear solver prior to calling SNESSolve().  In particular,
130c4421ceaSFande Kong      to employ an initial guess of zero, the user should explicitly set
131c4421ceaSFande Kong      this vector to zero by calling VecSet().
132c4421ceaSFande Kong   */
1339566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1349566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1359566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
136c4421ceaSFande Kong 
137c4421ceaSFande Kong   /*
138c4421ceaSFande Kong      Free work space.  All PETSc objects should be destroyed when they
139c4421ceaSFande Kong      are no longer needed.
140c4421ceaSFande Kong   */
1419371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1429371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1439371c9d4SSatish Balay   PetscCall(VecDestroy(&U));
1449371c9d4SSatish Balay   PetscCall(VecDestroy(&F));
1459371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1469371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1479566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
148b122ec5aSJacob Faibussowitsch   return 0;
149c4421ceaSFande Kong }
150f6dfbefdSBarry Smith 
151c4421ceaSFande Kong /*
152c4421ceaSFande Kong    FormInitialGuess - Computes initial guess.
153c4421ceaSFande Kong 
154c4421ceaSFande Kong    Input/Output Parameter:
155c4421ceaSFande Kong .  x - the solution vector
156c4421ceaSFande Kong */
FormInitialGuess(Vec x)157d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
158d71ae5a4SJacob Faibussowitsch {
159c4421ceaSFande Kong   PetscScalar pfive = .50;
1603ba16761SJacob Faibussowitsch 
1613ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1629566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
164c4421ceaSFande Kong }
165f6dfbefdSBarry Smith 
166c4421ceaSFande Kong /*
167c4421ceaSFande Kong    FormFunction - Evaluates nonlinear function, F(x).
168c4421ceaSFande Kong 
169c4421ceaSFande Kong    Input Parameters:
170c4421ceaSFande Kong .  snes - the SNES context
171c4421ceaSFande Kong .  x - input vector
172c4421ceaSFande Kong .  ctx - optional user-defined context, as set by SNESSetFunction()
173c4421ceaSFande Kong 
174c4421ceaSFande Kong    Output Parameter:
175c4421ceaSFande Kong .  f - function vector
176c4421ceaSFande Kong 
177c4421ceaSFande Kong    Note:
178c4421ceaSFande Kong    The user-defined context can contain any application-specific data
179c4421ceaSFande Kong    needed for the function evaluation (such as various parameters, work
180c4421ceaSFande Kong    vectors, and grid information).  In this program the context is just
181dd8e379bSPierre Jolivet    a vector containing the right-hand side of the discretized PDE.
182c4421ceaSFande Kong  */
183c4421ceaSFande Kong 
FormFunction(SNES snes,Vec x,Vec f,PetscCtx ctx)184*2a8381b2SBarry Smith PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, PetscCtx ctx)
185d71ae5a4SJacob Faibussowitsch {
186c4421ceaSFande Kong   Vec                g = (Vec)ctx;
187c4421ceaSFande Kong   const PetscScalar *xx, *gg;
188c4421ceaSFande Kong   PetscScalar       *ff, d;
189c4421ceaSFande Kong   PetscInt           i, n;
190c4421ceaSFande Kong 
1913ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
192c4421ceaSFande Kong   /*
193c4421ceaSFande Kong      Get pointers to vector data.
194c4421ceaSFande Kong        - For default PETSc vectors, VecGetArray() returns a pointer to
195c4421ceaSFande Kong          the data array.  Otherwise, the routine is implementation dependent.
196c4421ceaSFande Kong        - You MUST call VecRestoreArray() when you no longer need access to
197c4421ceaSFande Kong          the array.
198c4421ceaSFande Kong   */
1999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g, &gg));
202c4421ceaSFande Kong 
203c4421ceaSFande Kong   /*
204c4421ceaSFande Kong      Compute function
205c4421ceaSFande Kong   */
2069566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2079371c9d4SSatish Balay   d     = (PetscReal)(n - 1);
2089371c9d4SSatish Balay   d     = d * d;
209c4421ceaSFande Kong   ff[0] = xx[0];
210c4421ceaSFande 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];
211c4421ceaSFande Kong   ff[n - 1] = xx[n - 1] - 1.0;
212c4421ceaSFande Kong 
213c4421ceaSFande Kong   /*
214c4421ceaSFande Kong      Restore vectors
215c4421ceaSFande Kong   */
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
2189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g, &gg));
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
220c4421ceaSFande Kong }
221c4421ceaSFande Kong /* ------------------------------------------------------------------- */
222c4421ceaSFande Kong /*
223c4421ceaSFande Kong    FormJacobian - Evaluates Jacobian matrix.
224c4421ceaSFande Kong 
225c4421ceaSFande Kong    Input Parameters:
226c4421ceaSFande Kong .  snes - the SNES context
227c4421ceaSFande Kong .  x - input vector
228c4421ceaSFande Kong .  dummy - optional user-defined context (not used here)
229c4421ceaSFande Kong 
230c4421ceaSFande Kong    Output Parameters:
231c4421ceaSFande Kong .  jac - Jacobian matrix
2327addb90fSBarry Smith .  B - optionally different matrix used to construct the preconditioner
233c4421ceaSFande Kong 
234c4421ceaSFande Kong */
235c4421ceaSFande Kong 
FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void * dummy)236d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
237d71ae5a4SJacob Faibussowitsch {
238c4421ceaSFande Kong   const PetscScalar *xx;
239c4421ceaSFande Kong   PetscScalar        A[3], d;
240c4421ceaSFande Kong   PetscInt           i, n, j[3];
241c4421ceaSFande Kong 
2423ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
243c4421ceaSFande Kong   /*
244c4421ceaSFande Kong      Get pointer to vector data
245c4421ceaSFande Kong   */
2469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
247c4421ceaSFande Kong 
248c4421ceaSFande Kong   /*
249c4421ceaSFande Kong      Compute Jacobian entries and insert into matrix.
250c4421ceaSFande Kong       - Note that in this case we set all elements for a particular
251c4421ceaSFande Kong         row at once.
252c4421ceaSFande Kong   */
2539566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2549371c9d4SSatish Balay   d = (PetscReal)(n - 1);
2559371c9d4SSatish Balay   d = d * d;
256c4421ceaSFande Kong 
257c4421ceaSFande Kong   /*
258c4421ceaSFande Kong      Interior grid points
259c4421ceaSFande Kong   */
260c4421ceaSFande Kong   for (i = 1; i < n - 1; i++) {
2619371c9d4SSatish Balay     j[0] = i - 1;
2629371c9d4SSatish Balay     j[1] = i;
2639371c9d4SSatish Balay     j[2] = i + 1;
2649371c9d4SSatish Balay     A[0] = A[2] = d;
2659371c9d4SSatish Balay     A[1]        = -2.0 * d + 2.0 * xx[i];
2669566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
267c4421ceaSFande Kong   }
268c4421ceaSFande Kong 
269c4421ceaSFande Kong   /*
270c4421ceaSFande Kong      Boundary points
271c4421ceaSFande Kong   */
2729371c9d4SSatish Balay   i    = 0;
2739371c9d4SSatish Balay   A[0] = 1.0;
274c4421ceaSFande Kong 
2759566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
276c4421ceaSFande Kong 
2779371c9d4SSatish Balay   i    = n - 1;
2789371c9d4SSatish Balay   A[0] = 1.0;
279c4421ceaSFande Kong 
2809566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
281c4421ceaSFande Kong 
282c4421ceaSFande Kong   /*
283c4421ceaSFande Kong      Restore vector
284c4421ceaSFande Kong   */
2859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
286c4421ceaSFande Kong 
287c4421ceaSFande Kong   /*
288c4421ceaSFande Kong      Assemble matrix
289c4421ceaSFande Kong   */
2909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
292c4421ceaSFande Kong   if (jac != B) {
2939566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
295c4421ceaSFande Kong   }
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297c4421ceaSFande Kong }
298c4421ceaSFande Kong 
MySNESConvergedReasonView(SNES snes,PetscCtx ctx)299*2a8381b2SBarry Smith PetscErrorCode MySNESConvergedReasonView(SNES snes, PetscCtx ctx)
300d71ae5a4SJacob Faibussowitsch {
301c4421ceaSFande Kong   ReasonViewCtx      *monP   = (ReasonViewCtx *)ctx;
302c4421ceaSFande Kong   PetscViewer         viewer = monP->viewer;
303c4421ceaSFande Kong   SNESConvergedReason reason;
304c4421ceaSFande Kong   const char         *strreason;
305c4421ceaSFande Kong 
3063ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
3079566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes, &reason));
3089566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReasonString(snes, &strreason));
3099566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Customized SNES converged reason view\n"));
3109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, 1));
311c4421ceaSFande Kong   if (reason > 0) {
3129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Converged due to %s\n", strreason));
313c4421ceaSFande Kong   } else if (reason <= 0) {
3149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Did not converge due to %s\n", strreason));
315c4421ceaSFande Kong   }
3169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, 1));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318c4421ceaSFande Kong }
319c4421ceaSFande Kong 
MyKSPConvergedReasonView(KSP ksp,PetscCtx ctx)320*2a8381b2SBarry Smith PetscErrorCode MyKSPConvergedReasonView(KSP ksp, PetscCtx ctx)
321d71ae5a4SJacob Faibussowitsch {
322c4421ceaSFande Kong   ReasonViewCtx     *monP   = (ReasonViewCtx *)ctx;
323c4421ceaSFande Kong   PetscViewer        viewer = monP->viewer;
324c4421ceaSFande Kong   KSPConvergedReason reason;
325c4421ceaSFande Kong   const char        *reasonstr;
326c4421ceaSFande Kong 
3273ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
3289566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(ksp, &reason));
3299566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReasonString(ksp, &reasonstr));
3309566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, 2));
3319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Customized KSP converged reason view\n"));
3329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, 1));
333c4421ceaSFande Kong   if (reason > 0) {
3349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Converged due to %s\n", reasonstr));
335c4421ceaSFande Kong   } else if (reason <= 0) {
3369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Did not converge due to %s\n", reasonstr));
337c4421ceaSFande Kong   }
3389566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, 3));
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
340c4421ceaSFande Kong }
341c4421ceaSFande Kong 
342c4421ceaSFande Kong /*TEST
343c4421ceaSFande Kong 
344c4421ceaSFande Kong    test:
345c4421ceaSFande Kong       suffix: 1
346c4421ceaSFande Kong       nsize: 1
34736d43d94SBarry Smith       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
348c4421ceaSFande Kong 
349c4421ceaSFande Kong    test:
350c4421ceaSFande Kong       suffix: 2
351c4421ceaSFande Kong       nsize: 1
352c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel
35336d43d94SBarry Smith       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
354c4421ceaSFande Kong 
355c4421ceaSFande Kong    test:
356c4421ceaSFande Kong       suffix: 3
357c4421ceaSFande Kong       nsize: 1
358c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel -ksp_converged_reason
35936d43d94SBarry Smith       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
360c4421ceaSFande Kong 
361c4421ceaSFande Kong    test:
362c4421ceaSFande Kong       suffix: 4
363c4421ceaSFande Kong       nsize: 1
364c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel
36536d43d94SBarry Smith       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
366c4421ceaSFande Kong 
367c4421ceaSFande Kong    test:
368c4421ceaSFande Kong       suffix: 5
369c4421ceaSFande Kong       nsize: 1
370c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel -snes_converged_reason
37136d43d94SBarry Smith       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
372c4421ceaSFande Kong 
373c4421ceaSFande Kong TEST*/
374