xref: /petsc/src/snes/tutorials/ex6.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscInt       its,n = 5,i;
37c4421ceaSFande Kong   PetscMPIInt    size;
38c4421ceaSFande Kong   PetscScalar    h,xp,v;
39c4421ceaSFande Kong   MPI_Comm       comm;
40c4421ceaSFande Kong 
41*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
425f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
432c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!");
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
45c4421ceaSFande Kong   h    = 1.0/(n-1);
46c4421ceaSFande Kong   comm = PETSC_COMM_WORLD;
47c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48c4421ceaSFande Kong      Create nonlinear solver context
49c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50c4421ceaSFande Kong 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(comm,&snes));
52c4421ceaSFande Kong 
53c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54c4421ceaSFande Kong      Create vector data structures; set function evaluation routine
55c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56c4421ceaSFande Kong   /*
57c4421ceaSFande Kong      Note that we form 1 vector from scratch and then duplicate as needed.
58c4421ceaSFande Kong   */
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&x));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,n));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
65c4421ceaSFande Kong 
66c4421ceaSFande Kong   /*
67c4421ceaSFande Kong      Set function evaluation routine and vector
68c4421ceaSFande Kong   */
695f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormFunction,(void*)F));
70c4421ceaSFande Kong 
71c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72c4421ceaSFande Kong      Create matrix data structure; set Jacobian evaluation routine
73c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74c4421ceaSFande Kong 
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&J));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
79c4421ceaSFande Kong 
80c4421ceaSFande Kong   /*
81c4421ceaSFande Kong      Set Jacobian matrix data structure and default Jacobian evaluation
82c4421ceaSFande Kong      routine. User can override with:
83c4421ceaSFande Kong      -snes_fd : default finite differencing approximation of Jacobian
84c4421ceaSFande Kong      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
85c4421ceaSFande Kong                 (unless user explicitly sets preconditioner)
86c4421ceaSFande Kong      -snes_mf_operator : form preconditioning matrix as set by the user,
87c4421ceaSFande Kong                          but use matrix-free approx for Jacobian-vector
88c4421ceaSFande Kong                          products within Newton-Krylov method
89c4421ceaSFande Kong   */
90c4421ceaSFande Kong 
915f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
92c4421ceaSFande Kong 
93c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4421ceaSFande Kong      Customize nonlinear solver; set runtime options
95c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96c4421ceaSFande Kong 
97c4421ceaSFande Kong   /*
98c4421ceaSFande Kong      Set an optional user-defined reasonview routine
99c4421ceaSFande Kong   */
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIGetStdout(comm,&monP.viewer));
101c4421ceaSFande Kong   /* Just make sure we can not repeat addding the same function
102c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
103c4421ceaSFande Kong    */
104c4421ceaSFande Kong   for (i=0; i<4; i++) {
1055f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESConvergedReasonViewSet(snes,MySNESConvergedReasonView,&monP,0));
106c4421ceaSFande Kong   }
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
108c4421ceaSFande Kong   /* Just make sure we can not repeat addding the same function
109c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
110c4421ceaSFande Kong    */
111c4421ceaSFande Kong   for (i=0; i<4; i++) {
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPConvergedReasonViewSet(ksp,MyKSPConvergedReasonView,&monP,0));
113c4421ceaSFande Kong   }
114c4421ceaSFande Kong   /*
115c4421ceaSFande Kong      Set SNES/KSP/KSP/PC runtime options, e.g.,
116c4421ceaSFande Kong          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
117c4421ceaSFande Kong   */
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
119c4421ceaSFande Kong 
120c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4421ceaSFande Kong      Initialize application:
122c4421ceaSFande Kong      Store right-hand-side of PDE and exact solution
123c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4421ceaSFande Kong 
125c4421ceaSFande Kong   xp = 0.0;
126c4421ceaSFande Kong   for (i=0; i<n; i++) {
127c4421ceaSFande Kong     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES));
129c4421ceaSFande Kong     v    = xp*xp*xp;
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(U,1,&i,&v,INSERT_VALUES));
131c4421ceaSFande Kong     xp  += h;
132c4421ceaSFande Kong   }
133c4421ceaSFande Kong 
134c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135c4421ceaSFande Kong      Evaluate initial guess; then solve nonlinear system
136c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137c4421ceaSFande Kong   /*
138c4421ceaSFande Kong      Note: The user should initialize the vector, x, with the initial guess
139c4421ceaSFande Kong      for the nonlinear solver prior to calling SNESSolve().  In particular,
140c4421ceaSFande Kong      to employ an initial guess of zero, the user should explicitly set
141c4421ceaSFande Kong      this vector to zero by calling VecSet().
142c4421ceaSFande Kong   */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
146c4421ceaSFande Kong 
147c4421ceaSFande Kong   /*
148c4421ceaSFande Kong      Free work space.  All PETSc objects should be destroyed when they
149c4421ceaSFande Kong      are no longer needed.
150c4421ceaSFande Kong   */
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));  CHKERRQ(VecDestroy(&r));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));  CHKERRQ(VecDestroy(&F));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));  CHKERRQ(SNESDestroy(&snes));
1545f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscViewerDestroy(&monP.viewer));*/
155*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
156*b122ec5aSJacob Faibussowitsch   return 0;
157c4421ceaSFande Kong }
158c4421ceaSFande Kong /* ------------------------------------------------------------------- */
159c4421ceaSFande Kong /*
160c4421ceaSFande Kong    FormInitialGuess - Computes initial guess.
161c4421ceaSFande Kong 
162c4421ceaSFande Kong    Input/Output Parameter:
163c4421ceaSFande Kong .  x - the solution vector
164c4421ceaSFande Kong */
165c4421ceaSFande Kong PetscErrorCode FormInitialGuess(Vec x)
166c4421ceaSFande Kong {
167c4421ceaSFande Kong   PetscScalar    pfive = .50;
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
169c4421ceaSFande Kong   return 0;
170c4421ceaSFande Kong }
171c4421ceaSFande Kong /* ------------------------------------------------------------------- */
172c4421ceaSFande Kong /*
173c4421ceaSFande Kong    FormFunction - Evaluates nonlinear function, F(x).
174c4421ceaSFande Kong 
175c4421ceaSFande Kong    Input Parameters:
176c4421ceaSFande Kong .  snes - the SNES context
177c4421ceaSFande Kong .  x - input vector
178c4421ceaSFande Kong .  ctx - optional user-defined context, as set by SNESSetFunction()
179c4421ceaSFande Kong 
180c4421ceaSFande Kong    Output Parameter:
181c4421ceaSFande Kong .  f - function vector
182c4421ceaSFande Kong 
183c4421ceaSFande Kong    Note:
184c4421ceaSFande Kong    The user-defined context can contain any application-specific data
185c4421ceaSFande Kong    needed for the function evaluation (such as various parameters, work
186c4421ceaSFande Kong    vectors, and grid information).  In this program the context is just
187c4421ceaSFande Kong    a vector containing the right-hand-side of the discretized PDE.
188c4421ceaSFande Kong  */
189c4421ceaSFande Kong 
190c4421ceaSFande Kong PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
191c4421ceaSFande Kong {
192c4421ceaSFande Kong   Vec               g = (Vec)ctx;
193c4421ceaSFande Kong   const PetscScalar *xx,*gg;
194c4421ceaSFande Kong   PetscScalar       *ff,d;
195c4421ceaSFande Kong   PetscInt          i,n;
196c4421ceaSFande Kong 
197c4421ceaSFande Kong   /*
198c4421ceaSFande Kong      Get pointers to vector data.
199c4421ceaSFande Kong        - For default PETSc vectors, VecGetArray() returns a pointer to
200c4421ceaSFande Kong          the data array.  Otherwise, the routine is implementation dependent.
201c4421ceaSFande Kong        - You MUST call VecRestoreArray() when you no longer need access to
202c4421ceaSFande Kong          the array.
203c4421ceaSFande Kong   */
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&ff));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(g,&gg));
207c4421ceaSFande Kong 
208c4421ceaSFande Kong   /*
209c4421ceaSFande Kong      Compute function
210c4421ceaSFande Kong   */
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
212c4421ceaSFande Kong   d     = (PetscReal)(n - 1); d = d*d;
213c4421ceaSFande Kong   ff[0] = xx[0];
214c4421ceaSFande 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];
215c4421ceaSFande Kong   ff[n-1] = xx[n-1] - 1.0;
216c4421ceaSFande Kong 
217c4421ceaSFande Kong   /*
218c4421ceaSFande Kong      Restore vectors
219c4421ceaSFande Kong   */
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&ff));
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(g,&gg));
223c4421ceaSFande Kong   return 0;
224c4421ceaSFande Kong }
225c4421ceaSFande Kong /* ------------------------------------------------------------------- */
226c4421ceaSFande Kong /*
227c4421ceaSFande Kong    FormJacobian - Evaluates Jacobian matrix.
228c4421ceaSFande Kong 
229c4421ceaSFande Kong    Input Parameters:
230c4421ceaSFande Kong .  snes - the SNES context
231c4421ceaSFande Kong .  x - input vector
232c4421ceaSFande Kong .  dummy - optional user-defined context (not used here)
233c4421ceaSFande Kong 
234c4421ceaSFande Kong    Output Parameters:
235c4421ceaSFande Kong .  jac - Jacobian matrix
236c4421ceaSFande Kong .  B - optionally different preconditioning matrix
237c4421ceaSFande Kong 
238c4421ceaSFande Kong */
239c4421ceaSFande Kong 
240c4421ceaSFande Kong PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
241c4421ceaSFande Kong {
242c4421ceaSFande Kong   const PetscScalar *xx;
243c4421ceaSFande Kong   PetscScalar       A[3],d;
244c4421ceaSFande Kong   PetscInt          i,n,j[3];
245c4421ceaSFande Kong 
246c4421ceaSFande Kong   /*
247c4421ceaSFande Kong      Get pointer to vector data
248c4421ceaSFande Kong   */
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
250c4421ceaSFande Kong 
251c4421ceaSFande Kong   /*
252c4421ceaSFande Kong      Compute Jacobian entries and insert into matrix.
253c4421ceaSFande Kong       - Note that in this case we set all elements for a particular
254c4421ceaSFande Kong         row at once.
255c4421ceaSFande Kong   */
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
257c4421ceaSFande Kong   d    = (PetscReal)(n - 1); d = d*d;
258c4421ceaSFande Kong 
259c4421ceaSFande Kong   /*
260c4421ceaSFande Kong      Interior grid points
261c4421ceaSFande Kong   */
262c4421ceaSFande Kong   for (i=1; i<n-1; i++) {
263c4421ceaSFande Kong     j[0] = i - 1; j[1] = i; j[2] = i + 1;
264c4421ceaSFande Kong     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
2655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
266c4421ceaSFande Kong   }
267c4421ceaSFande Kong 
268c4421ceaSFande Kong   /*
269c4421ceaSFande Kong      Boundary points
270c4421ceaSFande Kong   */
271c4421ceaSFande Kong   i = 0;   A[0] = 1.0;
272c4421ceaSFande Kong 
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
274c4421ceaSFande Kong 
275c4421ceaSFande Kong   i = n-1; A[0] = 1.0;
276c4421ceaSFande Kong 
2775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
278c4421ceaSFande Kong 
279c4421ceaSFande Kong   /*
280c4421ceaSFande Kong      Restore vector
281c4421ceaSFande Kong   */
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
283c4421ceaSFande Kong 
284c4421ceaSFande Kong   /*
285c4421ceaSFande Kong      Assemble matrix
286c4421ceaSFande Kong   */
2875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
289c4421ceaSFande Kong   if (jac != B) {
2905f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
2915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
292c4421ceaSFande Kong   }
293c4421ceaSFande Kong   return 0;
294c4421ceaSFande Kong }
295c4421ceaSFande Kong 
296c4421ceaSFande Kong PetscErrorCode MySNESConvergedReasonView(SNES snes,void *ctx)
297c4421ceaSFande Kong {
298c4421ceaSFande Kong   ReasonViewCtx         *monP = (ReasonViewCtx*) ctx;
299c4421ceaSFande Kong   PetscViewer           viewer = monP->viewer;
300c4421ceaSFande Kong   SNESConvergedReason   reason;
301c4421ceaSFande Kong   const char            *strreason;
302c4421ceaSFande Kong 
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetConvergedReason(snes,&reason));
3045f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetConvergedReasonString(snes,&strreason));
3055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Customized SNES converged reason view\n"));
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIAddTab(viewer,1));
307c4421ceaSFande Kong   if (reason > 0) {
3085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",strreason));
309c4421ceaSFande Kong   } else if (reason <= 0) {
3105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",strreason));
311c4421ceaSFande Kong   }
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIISubtractTab(viewer,1));
313c4421ceaSFande Kong   return 0;
314c4421ceaSFande Kong }
315c4421ceaSFande Kong 
316c4421ceaSFande Kong PetscErrorCode MyKSPConvergedReasonView(KSP ksp,void *ctx)
317c4421ceaSFande Kong {
318c4421ceaSFande Kong   ReasonViewCtx         *monP = (ReasonViewCtx*) ctx;
319c4421ceaSFande Kong   PetscViewer           viewer = monP->viewer;
320c4421ceaSFande Kong   KSPConvergedReason    reason;
321c4421ceaSFande Kong   const char            *reasonstr;
322c4421ceaSFande Kong 
3235f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetConvergedReason(ksp,&reason));
3245f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetConvergedReasonString(ksp,&reasonstr));
3255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIAddTab(viewer,2));
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Customized KSP converged reason view\n"));
3275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIAddTab(viewer,1));
328c4421ceaSFande Kong   if (reason > 0) {
3295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",reasonstr));
330c4421ceaSFande Kong   } else if (reason <= 0) {
3315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",reasonstr));
332c4421ceaSFande Kong   }
3335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIISubtractTab(viewer,3));
334c4421ceaSFande Kong   return 0;
335c4421ceaSFande Kong }
336c4421ceaSFande Kong 
337c4421ceaSFande Kong /*TEST
338c4421ceaSFande Kong 
339c4421ceaSFande Kong    test:
340c4421ceaSFande Kong       suffix: 1
341c4421ceaSFande Kong       nsize: 1
342bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
343c4421ceaSFande Kong 
344c4421ceaSFande Kong    test:
345c4421ceaSFande Kong       suffix: 2
346c4421ceaSFande Kong       nsize: 1
347c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel
348bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
349c4421ceaSFande Kong 
350c4421ceaSFande Kong    test:
351c4421ceaSFande Kong       suffix: 3
352c4421ceaSFande Kong       nsize: 1
353c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel -ksp_converged_reason
354bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
355c4421ceaSFande Kong 
356c4421ceaSFande Kong    test:
357c4421ceaSFande Kong       suffix: 4
358c4421ceaSFande Kong       nsize: 1
359c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel
360bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
361c4421ceaSFande Kong 
362c4421ceaSFande Kong    test:
363c4421ceaSFande Kong       suffix: 5
364c4421ceaSFande Kong       nsize: 1
365c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel -snes_converged_reason
366bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
367c4421ceaSFande Kong 
368c4421ceaSFande Kong TEST*/
369