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