xref: /petsc/src/snes/tutorials/ex6.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 /*T
6c4421ceaSFande Kong    Concepts: SNES^basic uniprocessor example
7c4421ceaSFande Kong    Concepts: SNES^setting a user-defined reasonview routine
8c4421ceaSFande Kong    Processors: 1
9c4421ceaSFande Kong T*/
10c4421ceaSFande Kong 
11c4421ceaSFande Kong #include <petscsnes.h>
12c4421ceaSFande Kong 
13c4421ceaSFande Kong /*
14c4421ceaSFande Kong    User-defined routines
15c4421ceaSFande Kong */
16c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
17c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
18c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormInitialGuess(Vec);
19c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MySNESConvergedReasonView(SNES,void*);
20c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MyKSPConvergedReasonView(KSP,void*);
21c4421ceaSFande Kong 
22c4421ceaSFande Kong /*
23c4421ceaSFande Kong    User-defined context for monitoring
24c4421ceaSFande Kong */
25c4421ceaSFande Kong typedef struct {
26c4421ceaSFande Kong   PetscViewer viewer;
27c4421ceaSFande Kong } ReasonViewCtx;
28c4421ceaSFande Kong 
29c4421ceaSFande Kong int main(int argc,char **argv)
30c4421ceaSFande Kong {
31c4421ceaSFande Kong   SNES           snes;                /* SNES context */
32c4421ceaSFande Kong   KSP            ksp;                 /* KSP context */
33c4421ceaSFande Kong   Vec            x,r,F,U;             /* vectors */
34c4421ceaSFande Kong   Mat            J;                   /* Jacobian matrix */
35c4421ceaSFande Kong   ReasonViewCtx     monP;             /* monitoring context */
36c4421ceaSFande Kong   PetscErrorCode ierr;
37c4421ceaSFande Kong   PetscInt       its,n = 5,i;
38c4421ceaSFande Kong   PetscMPIInt    size;
39c4421ceaSFande Kong   PetscScalar    h,xp,v;
40c4421ceaSFande Kong   MPI_Comm       comm;
41c4421ceaSFande Kong 
42c4421ceaSFande Kong   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
43*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
442c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!");
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
46c4421ceaSFande Kong   h    = 1.0/(n-1);
47c4421ceaSFande Kong   comm = PETSC_COMM_WORLD;
48c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49c4421ceaSFande Kong      Create nonlinear solver context
50c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51c4421ceaSFande Kong 
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(comm,&snes));
53c4421ceaSFande Kong 
54c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55c4421ceaSFande Kong      Create vector data structures; set function evaluation routine
56c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57c4421ceaSFande Kong   /*
58c4421ceaSFande Kong      Note that we form 1 vector from scratch and then duplicate as needed.
59c4421ceaSFande Kong   */
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&x));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,n));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
66c4421ceaSFande Kong 
67c4421ceaSFande Kong   /*
68c4421ceaSFande Kong      Set function evaluation routine and vector
69c4421ceaSFande Kong   */
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormFunction,(void*)F));
71c4421ceaSFande Kong 
72c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73c4421ceaSFande Kong      Create matrix data structure; set Jacobian evaluation routine
74c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75c4421ceaSFande Kong 
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&J));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
80c4421ceaSFande Kong 
81c4421ceaSFande Kong   /*
82c4421ceaSFande Kong      Set Jacobian matrix data structure and default Jacobian evaluation
83c4421ceaSFande Kong      routine. User can override with:
84c4421ceaSFande Kong      -snes_fd : default finite differencing approximation of Jacobian
85c4421ceaSFande Kong      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
86c4421ceaSFande Kong                 (unless user explicitly sets preconditioner)
87c4421ceaSFande Kong      -snes_mf_operator : form preconditioning matrix as set by the user,
88c4421ceaSFande Kong                          but use matrix-free approx for Jacobian-vector
89c4421ceaSFande Kong                          products within Newton-Krylov method
90c4421ceaSFande Kong   */
91c4421ceaSFande Kong 
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
93c4421ceaSFande Kong 
94c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4421ceaSFande Kong      Customize nonlinear solver; set runtime options
96c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97c4421ceaSFande Kong 
98c4421ceaSFande Kong   /*
99c4421ceaSFande Kong      Set an optional user-defined reasonview routine
100c4421ceaSFande Kong   */
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIGetStdout(comm,&monP.viewer));
102c4421ceaSFande Kong   /* Just make sure we can not repeat addding the same function
103c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
104c4421ceaSFande Kong    */
105c4421ceaSFande Kong   for (i=0; i<4; i++) {
106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESConvergedReasonViewSet(snes,MySNESConvergedReasonView,&monP,0));
107c4421ceaSFande Kong   }
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
109c4421ceaSFande Kong   /* Just make sure we can not repeat addding the same function
110c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
111c4421ceaSFande Kong    */
112c4421ceaSFande Kong   for (i=0; i<4; i++) {
113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPConvergedReasonViewSet(ksp,MyKSPConvergedReasonView,&monP,0));
114c4421ceaSFande Kong   }
115c4421ceaSFande Kong   /*
116c4421ceaSFande Kong      Set SNES/KSP/KSP/PC runtime options, e.g.,
117c4421ceaSFande Kong          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
118c4421ceaSFande Kong   */
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
120c4421ceaSFande Kong 
121c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4421ceaSFande Kong      Initialize application:
123c4421ceaSFande Kong      Store right-hand-side of PDE and exact solution
124c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125c4421ceaSFande Kong 
126c4421ceaSFande Kong   xp = 0.0;
127c4421ceaSFande Kong   for (i=0; i<n; i++) {
128c4421ceaSFande Kong     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES));
130c4421ceaSFande Kong     v    = xp*xp*xp;
131*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(U,1,&i,&v,INSERT_VALUES));
132c4421ceaSFande Kong     xp  += h;
133c4421ceaSFande Kong   }
134c4421ceaSFande Kong 
135c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4421ceaSFande Kong      Evaluate initial guess; then solve nonlinear system
137c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138c4421ceaSFande Kong   /*
139c4421ceaSFande Kong      Note: The user should initialize the vector, x, with the initial guess
140c4421ceaSFande Kong      for the nonlinear solver prior to calling SNESSolve().  In particular,
141c4421ceaSFande Kong      to employ an initial guess of zero, the user should explicitly set
142c4421ceaSFande Kong      this vector to zero by calling VecSet().
143c4421ceaSFande Kong   */
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
147c4421ceaSFande Kong 
148c4421ceaSFande Kong   /*
149c4421ceaSFande Kong      Free work space.  All PETSc objects should be destroyed when they
150c4421ceaSFande Kong      are no longer needed.
151c4421ceaSFande Kong   */
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));  CHKERRQ(VecDestroy(&r));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));  CHKERRQ(VecDestroy(&F));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));  CHKERRQ(SNESDestroy(&snes));
155*5f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscViewerDestroy(&monP.viewer));*/
156c4421ceaSFande Kong   ierr = PetscFinalize();
157c4421ceaSFande Kong   return ierr;
158c4421ceaSFande Kong }
159c4421ceaSFande Kong /* ------------------------------------------------------------------- */
160c4421ceaSFande Kong /*
161c4421ceaSFande Kong    FormInitialGuess - Computes initial guess.
162c4421ceaSFande Kong 
163c4421ceaSFande Kong    Input/Output Parameter:
164c4421ceaSFande Kong .  x - the solution vector
165c4421ceaSFande Kong */
166c4421ceaSFande Kong PetscErrorCode FormInitialGuess(Vec x)
167c4421ceaSFande Kong {
168c4421ceaSFande Kong   PetscScalar    pfive = .50;
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
170c4421ceaSFande Kong   return 0;
171c4421ceaSFande Kong }
172c4421ceaSFande Kong /* ------------------------------------------------------------------- */
173c4421ceaSFande Kong /*
174c4421ceaSFande Kong    FormFunction - Evaluates nonlinear function, F(x).
175c4421ceaSFande Kong 
176c4421ceaSFande Kong    Input Parameters:
177c4421ceaSFande Kong .  snes - the SNES context
178c4421ceaSFande Kong .  x - input vector
179c4421ceaSFande Kong .  ctx - optional user-defined context, as set by SNESSetFunction()
180c4421ceaSFande Kong 
181c4421ceaSFande Kong    Output Parameter:
182c4421ceaSFande Kong .  f - function vector
183c4421ceaSFande Kong 
184c4421ceaSFande Kong    Note:
185c4421ceaSFande Kong    The user-defined context can contain any application-specific data
186c4421ceaSFande Kong    needed for the function evaluation (such as various parameters, work
187c4421ceaSFande Kong    vectors, and grid information).  In this program the context is just
188c4421ceaSFande Kong    a vector containing the right-hand-side of the discretized PDE.
189c4421ceaSFande Kong  */
190c4421ceaSFande Kong 
191c4421ceaSFande Kong PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
192c4421ceaSFande Kong {
193c4421ceaSFande Kong   Vec               g = (Vec)ctx;
194c4421ceaSFande Kong   const PetscScalar *xx,*gg;
195c4421ceaSFande Kong   PetscScalar       *ff,d;
196c4421ceaSFande Kong   PetscInt          i,n;
197c4421ceaSFande Kong 
198c4421ceaSFande Kong   /*
199c4421ceaSFande Kong      Get pointers to vector data.
200c4421ceaSFande Kong        - For default PETSc vectors, VecGetArray() returns a pointer to
201c4421ceaSFande Kong          the data array.  Otherwise, the routine is implementation dependent.
202c4421ceaSFande Kong        - You MUST call VecRestoreArray() when you no longer need access to
203c4421ceaSFande Kong          the array.
204c4421ceaSFande Kong   */
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&ff));
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(g,&gg));
208c4421ceaSFande Kong 
209c4421ceaSFande Kong   /*
210c4421ceaSFande Kong      Compute function
211c4421ceaSFande Kong   */
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
213c4421ceaSFande Kong   d     = (PetscReal)(n - 1); d = d*d;
214c4421ceaSFande Kong   ff[0] = xx[0];
215c4421ceaSFande 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];
216c4421ceaSFande Kong   ff[n-1] = xx[n-1] - 1.0;
217c4421ceaSFande Kong 
218c4421ceaSFande Kong   /*
219c4421ceaSFande Kong      Restore vectors
220c4421ceaSFande Kong   */
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&ff));
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(g,&gg));
224c4421ceaSFande Kong   return 0;
225c4421ceaSFande Kong }
226c4421ceaSFande Kong /* ------------------------------------------------------------------- */
227c4421ceaSFande Kong /*
228c4421ceaSFande Kong    FormJacobian - Evaluates Jacobian matrix.
229c4421ceaSFande Kong 
230c4421ceaSFande Kong    Input Parameters:
231c4421ceaSFande Kong .  snes - the SNES context
232c4421ceaSFande Kong .  x - input vector
233c4421ceaSFande Kong .  dummy - optional user-defined context (not used here)
234c4421ceaSFande Kong 
235c4421ceaSFande Kong    Output Parameters:
236c4421ceaSFande Kong .  jac - Jacobian matrix
237c4421ceaSFande Kong .  B - optionally different preconditioning matrix
238c4421ceaSFande Kong 
239c4421ceaSFande Kong */
240c4421ceaSFande Kong 
241c4421ceaSFande Kong PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
242c4421ceaSFande Kong {
243c4421ceaSFande Kong   const PetscScalar *xx;
244c4421ceaSFande Kong   PetscScalar       A[3],d;
245c4421ceaSFande Kong   PetscInt          i,n,j[3];
246c4421ceaSFande Kong 
247c4421ceaSFande Kong   /*
248c4421ceaSFande Kong      Get pointer to vector data
249c4421ceaSFande Kong   */
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
251c4421ceaSFande Kong 
252c4421ceaSFande Kong   /*
253c4421ceaSFande Kong      Compute Jacobian entries and insert into matrix.
254c4421ceaSFande Kong       - Note that in this case we set all elements for a particular
255c4421ceaSFande Kong         row at once.
256c4421ceaSFande Kong   */
257*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
258c4421ceaSFande Kong   d    = (PetscReal)(n - 1); d = d*d;
259c4421ceaSFande Kong 
260c4421ceaSFande Kong   /*
261c4421ceaSFande Kong      Interior grid points
262c4421ceaSFande Kong   */
263c4421ceaSFande Kong   for (i=1; i<n-1; i++) {
264c4421ceaSFande Kong     j[0] = i - 1; j[1] = i; j[2] = i + 1;
265c4421ceaSFande Kong     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
266*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
267c4421ceaSFande Kong   }
268c4421ceaSFande Kong 
269c4421ceaSFande Kong   /*
270c4421ceaSFande Kong      Boundary points
271c4421ceaSFande Kong   */
272c4421ceaSFande Kong   i = 0;   A[0] = 1.0;
273c4421ceaSFande Kong 
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
275c4421ceaSFande Kong 
276c4421ceaSFande Kong   i = n-1; A[0] = 1.0;
277c4421ceaSFande Kong 
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
279c4421ceaSFande Kong 
280c4421ceaSFande Kong   /*
281c4421ceaSFande Kong      Restore vector
282c4421ceaSFande Kong   */
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
284c4421ceaSFande Kong 
285c4421ceaSFande Kong   /*
286c4421ceaSFande Kong      Assemble matrix
287c4421ceaSFande Kong   */
288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
290c4421ceaSFande Kong   if (jac != B) {
291*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
292*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
293c4421ceaSFande Kong   }
294c4421ceaSFande Kong   return 0;
295c4421ceaSFande Kong }
296c4421ceaSFande Kong 
297c4421ceaSFande Kong PetscErrorCode MySNESConvergedReasonView(SNES snes,void *ctx)
298c4421ceaSFande Kong {
299c4421ceaSFande Kong   ReasonViewCtx         *monP = (ReasonViewCtx*) ctx;
300c4421ceaSFande Kong   PetscViewer           viewer = monP->viewer;
301c4421ceaSFande Kong   SNESConvergedReason   reason;
302c4421ceaSFande Kong   const char            *strreason;
303c4421ceaSFande Kong 
304*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetConvergedReason(snes,&reason));
305*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetConvergedReasonString(snes,&strreason));
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Customized SNES converged reason view\n"));
307*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIAddTab(viewer,1));
308c4421ceaSFande Kong   if (reason > 0) {
309*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",strreason));
310c4421ceaSFande Kong   } else if (reason <= 0) {
311*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",strreason));
312c4421ceaSFande Kong   }
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIISubtractTab(viewer,1));
314c4421ceaSFande Kong   return 0;
315c4421ceaSFande Kong }
316c4421ceaSFande Kong 
317c4421ceaSFande Kong PetscErrorCode MyKSPConvergedReasonView(KSP ksp,void *ctx)
318c4421ceaSFande Kong {
319c4421ceaSFande Kong   ReasonViewCtx         *monP = (ReasonViewCtx*) ctx;
320c4421ceaSFande Kong   PetscViewer           viewer = monP->viewer;
321c4421ceaSFande Kong   KSPConvergedReason    reason;
322c4421ceaSFande Kong   const char            *reasonstr;
323c4421ceaSFande Kong 
324*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetConvergedReason(ksp,&reason));
325*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetConvergedReasonString(ksp,&reasonstr));
326*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIAddTab(viewer,2));
327*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Customized KSP converged reason view\n"));
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIAddTab(viewer,1));
329c4421ceaSFande Kong   if (reason > 0) {
330*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",reasonstr));
331c4421ceaSFande Kong   } else if (reason <= 0) {
332*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",reasonstr));
333c4421ceaSFande Kong   }
334*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIISubtractTab(viewer,3));
335c4421ceaSFande Kong   return 0;
336c4421ceaSFande Kong }
337c4421ceaSFande Kong 
338c4421ceaSFande Kong /*TEST
339c4421ceaSFande Kong 
340c4421ceaSFande Kong    test:
341c4421ceaSFande Kong       suffix: 1
342c4421ceaSFande Kong       nsize: 1
343bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
344c4421ceaSFande Kong 
345c4421ceaSFande Kong    test:
346c4421ceaSFande Kong       suffix: 2
347c4421ceaSFande Kong       nsize: 1
348c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel
349bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
350c4421ceaSFande Kong 
351c4421ceaSFande Kong    test:
352c4421ceaSFande Kong       suffix: 3
353c4421ceaSFande Kong       nsize: 1
354c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel -ksp_converged_reason
355bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
356c4421ceaSFande Kong 
357c4421ceaSFande Kong    test:
358c4421ceaSFande Kong       suffix: 4
359c4421ceaSFande Kong       nsize: 1
360c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel
361bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
362c4421ceaSFande Kong 
363c4421ceaSFande Kong    test:
364c4421ceaSFande Kong       suffix: 5
365c4421ceaSFande Kong       nsize: 1
366c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel -snes_converged_reason
367bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
368c4421ceaSFande Kong 
369c4421ceaSFande Kong TEST*/
370