xref: /petsc/src/snes/tutorials/ex3.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1 static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel.\n\
2 This example employs a user-defined monitoring routine and optionally a user-defined\n\
3 routine to check candidate iterates produced by line search routines.\n\
4 The command line options include:\n\
5   -pre_check_iterates : activate checking of iterates\n\
6   -post_check_iterates : activate checking of iterates\n\
7   -check_tol <tol>: set tolerance for iterate checking\n\
8   -user_precond : activate a (trivial) user-defined preconditioner\n\n";
9 
10 /*T
11    Concepts: SNES^basic parallel example
12    Concepts: SNES^setting a user-defined monitoring routine
13    Processors: n
14 T*/
15 
16 /*
17    Include "petscdm.h" so that we can use data management objects (DMs)
18    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
19    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
20    file automatically includes:
21      petscsys.h    - base PETSc routines
22      petscvec.h    - vectors
23      petscmat.h    - matrices
24      petscis.h     - index sets
25      petscksp.h    - Krylov subspace methods
26      petscviewer.h - viewers
27      petscpc.h     - preconditioners
28      petscksp.h    - linear solvers
29 */
30 
31 #include <petscdm.h>
32 #include <petscdmda.h>
33 #include <petscsnes.h>
34 
35 /*
36    User-defined routines.
37 */
38 PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
39 PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
40 PetscErrorCode FormInitialGuess(Vec);
41 PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
42 PetscErrorCode PreCheck(SNESLineSearch,Vec,Vec,PetscBool*,void*);
43 PetscErrorCode PostCheck(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
44 PetscErrorCode PostSetSubKSP(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
45 PetscErrorCode MatrixFreePreconditioner(PC,Vec,Vec);
46 
47 /*
48    User-defined application context
49 */
50 typedef struct {
51   DM          da;      /* distributed array */
52   Vec         F;       /* right-hand-side of PDE */
53   PetscMPIInt rank;    /* rank of processor */
54   PetscMPIInt size;    /* size of communicator */
55   PetscReal   h;       /* mesh spacing */
56   PetscBool   sjerr;   /* if or not to test jacobian domain error */
57 } ApplicationCtx;
58 
59 /*
60    User-defined context for monitoring
61 */
62 typedef struct {
63   PetscViewer viewer;
64 } MonitorCtx;
65 
66 /*
67    User-defined context for checking candidate iterates that are
68    determined by line search methods
69 */
70 typedef struct {
71   Vec            last_step;  /* previous iterate */
72   PetscReal      tolerance;  /* tolerance for changes between successive iterates */
73   ApplicationCtx *user;
74 } StepCheckCtx;
75 
76 typedef struct {
77   PetscInt its0; /* num of prevous outer KSP iterations */
78 } SetSubKSPCtx;
79 
80 int main(int argc,char **argv)
81 {
82   SNES           snes;                 /* SNES context */
83   SNESLineSearch linesearch;           /* SNESLineSearch context */
84   Mat            J;                    /* Jacobian matrix */
85   ApplicationCtx ctx;                  /* user-defined context */
86   Vec            x,r,U,F;              /* vectors */
87   MonitorCtx     monP;                 /* monitoring context */
88   StepCheckCtx   checkP;               /* step-checking context */
89   SetSubKSPCtx   checkP1;
90   PetscBool      pre_check,post_check,post_setsubksp; /* flag indicating whether we're checking candidate iterates */
91   PetscScalar    xp,*FF,*UU,none = -1.0;
92   PetscErrorCode ierr;
93   PetscInt       its,N = 5,i,maxit,maxf,xs,xm;
94   PetscReal      abstol,rtol,stol,norm;
95   PetscBool      flg,viewinitial = PETSC_FALSE;
96 
97   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
98   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank));
99   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size));
100   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
101   ctx.h = 1.0/(N-1);
102   ctx.sjerr = PETSC_FALSE;
103   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_jacobian_domain_error",&ctx.sjerr,NULL));
104   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL));
105 
106   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107      Create nonlinear solver context
108      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109 
110   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
111 
112   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113      Create vector data structures; set function evaluation routine
114      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115 
116   /*
117      Create distributed array (DMDA) to manage parallel grid and vectors
118   */
119   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da));
120   CHKERRQ(DMSetFromOptions(ctx.da));
121   CHKERRQ(DMSetUp(ctx.da));
122 
123   /*
124      Extract global and local vectors from DMDA; then duplicate for remaining
125      vectors that are the same types
126   */
127   CHKERRQ(DMCreateGlobalVector(ctx.da,&x));
128   CHKERRQ(VecDuplicate(x,&r));
129   CHKERRQ(VecDuplicate(x,&F)); ctx.F = F;
130   CHKERRQ(VecDuplicate(x,&U));
131 
132   /*
133      Set function evaluation routine and vector.  Whenever the nonlinear
134      solver needs to compute the nonlinear function, it will call this
135      routine.
136       - Note that the final routine argument is the user-defined
137         context that provides application-specific data for the
138         function evaluation routine.
139   */
140   CHKERRQ(SNESSetFunction(snes,r,FormFunction,&ctx));
141 
142   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143      Create matrix data structure; set Jacobian evaluation routine
144      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145 
146   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
147   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N));
148   CHKERRQ(MatSetFromOptions(J));
149   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
150   CHKERRQ(MatMPIAIJSetPreallocation(J,3,NULL,3,NULL));
151 
152   /*
153      Set Jacobian matrix data structure and default Jacobian evaluation
154      routine.  Whenever the nonlinear solver needs to compute the
155      Jacobian matrix, it will call this routine.
156       - Note that the final routine argument is the user-defined
157         context that provides application-specific data for the
158         Jacobian evaluation routine.
159   */
160   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,&ctx));
161 
162   /*
163      Optionally allow user-provided preconditioner
164    */
165   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-user_precond",&flg));
166   if (flg) {
167     KSP ksp;
168     PC  pc;
169     CHKERRQ(SNESGetKSP(snes,&ksp));
170     CHKERRQ(KSPGetPC(ksp,&pc));
171     CHKERRQ(PCSetType(pc,PCSHELL));
172     CHKERRQ(PCShellSetApply(pc,MatrixFreePreconditioner));
173   }
174 
175   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176      Customize nonlinear solver; set runtime options
177    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178 
179   /*
180      Set an optional user-defined monitoring routine
181   */
182   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer));
183   CHKERRQ(SNESMonitorSet(snes,Monitor,&monP,0));
184 
185   /*
186      Set names for some vectors to facilitate monitoring (optional)
187   */
188   CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
189   CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution"));
190 
191   /*
192      Set SNES/KSP/KSP/PC runtime options, e.g.,
193          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
194   */
195   CHKERRQ(SNESSetFromOptions(snes));
196 
197   /*
198      Set an optional user-defined routine to check the validity of candidate
199      iterates that are determined by line search methods
200   */
201   CHKERRQ(SNESGetLineSearch(snes, &linesearch));
202   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-post_check_iterates",&post_check));
203 
204   if (post_check) {
205     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating post step checking routine\n"));
206     CHKERRQ(SNESLineSearchSetPostCheck(linesearch,PostCheck,&checkP));
207     CHKERRQ(VecDuplicate(x,&(checkP.last_step)));
208 
209     checkP.tolerance = 1.0;
210     checkP.user      = &ctx;
211 
212     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-check_tol",&checkP.tolerance,NULL));
213   }
214 
215   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-post_setsubksp",&post_setsubksp));
216   if (post_setsubksp) {
217     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating post setsubksp\n"));
218     CHKERRQ(SNESLineSearchSetPostCheck(linesearch,PostSetSubKSP,&checkP1));
219   }
220 
221   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-pre_check_iterates",&pre_check));
222   if (pre_check) {
223     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating pre step checking routine\n"));
224     CHKERRQ(SNESLineSearchSetPreCheck(linesearch,PreCheck,&checkP));
225   }
226 
227   /*
228      Print parameters used for convergence testing (optional) ... just
229      to demonstrate this routine; this information is also printed with
230      the option -snes_view
231   */
232   CHKERRQ(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
233   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
234 
235   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236      Initialize application:
237      Store right-hand-side of PDE and exact solution
238    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239 
240   /*
241      Get local grid boundaries (for 1-dimensional DMDA):
242        xs, xm - starting grid index, width of local grid (no ghost points)
243   */
244   CHKERRQ(DMDAGetCorners(ctx.da,&xs,NULL,NULL,&xm,NULL,NULL));
245 
246   /*
247      Get pointers to vector data
248   */
249   CHKERRQ(DMDAVecGetArray(ctx.da,F,&FF));
250   CHKERRQ(DMDAVecGetArray(ctx.da,U,&UU));
251 
252   /*
253      Compute local vector entries
254   */
255   xp = ctx.h*xs;
256   for (i=xs; i<xs+xm; i++) {
257     FF[i] = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
258     UU[i] = xp*xp*xp;
259     xp   += ctx.h;
260   }
261 
262   /*
263      Restore vectors
264   */
265   CHKERRQ(DMDAVecRestoreArray(ctx.da,F,&FF));
266   CHKERRQ(DMDAVecRestoreArray(ctx.da,U,&UU));
267   if (viewinitial) {
268     CHKERRQ(VecView(U,0));
269     CHKERRQ(VecView(F,0));
270   }
271 
272   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
273      Evaluate initial guess; then solve nonlinear system
274    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
275 
276   /*
277      Note: The user should initialize the vector, x, with the initial guess
278      for the nonlinear solver prior to calling SNESSolve().  In particular,
279      to employ an initial guess of zero, the user should explicitly set
280      this vector to zero by calling VecSet().
281   */
282   CHKERRQ(FormInitialGuess(x));
283   CHKERRQ(SNESSolve(snes,NULL,x));
284   CHKERRQ(SNESGetIterationNumber(snes,&its));
285   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
286 
287   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288      Check solution and clean up
289    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
290 
291   /*
292      Check the error
293   */
294   CHKERRQ(VecAXPY(x,none,U));
295   CHKERRQ(VecNorm(x,NORM_2,&norm));
296   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its));
297   if (ctx.sjerr) {
298     SNESType       snestype;
299     CHKERRQ(SNESGetType(snes,&snestype));
300     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"SNES Type %s\n",snestype));
301   }
302 
303   /*
304      Free work space.  All PETSc objects should be destroyed when they
305      are no longer needed.
306   */
307   CHKERRQ(PetscViewerDestroy(&monP.viewer));
308   if (post_check) CHKERRQ(VecDestroy(&checkP.last_step));
309   CHKERRQ(VecDestroy(&x));
310   CHKERRQ(VecDestroy(&r));
311   CHKERRQ(VecDestroy(&U));
312   CHKERRQ(VecDestroy(&F));
313   CHKERRQ(MatDestroy(&J));
314   CHKERRQ(SNESDestroy(&snes));
315   CHKERRQ(DMDestroy(&ctx.da));
316   ierr = PetscFinalize();
317   return ierr;
318 }
319 
320 /* ------------------------------------------------------------------- */
321 /*
322    FormInitialGuess - Computes initial guess.
323 
324    Input/Output Parameter:
325 .  x - the solution vector
326 */
327 PetscErrorCode FormInitialGuess(Vec x)
328 {
329   PetscScalar    pfive = .50;
330 
331   PetscFunctionBeginUser;
332   CHKERRQ(VecSet(x,pfive));
333   PetscFunctionReturn(0);
334 }
335 
336 /* ------------------------------------------------------------------- */
337 /*
338    FormFunction - Evaluates nonlinear function, F(x).
339 
340    Input Parameters:
341 .  snes - the SNES context
342 .  x - input vector
343 .  ctx - optional user-defined context, as set by SNESSetFunction()
344 
345    Output Parameter:
346 .  f - function vector
347 
348    Note:
349    The user-defined context can contain any application-specific
350    data needed for the function evaluation.
351 */
352 PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
353 {
354   ApplicationCtx *user = (ApplicationCtx*) ctx;
355   DM             da    = user->da;
356   PetscScalar    *xx,*ff,*FF,d;
357   PetscInt       i,M,xs,xm;
358   Vec            xlocal;
359 
360   PetscFunctionBeginUser;
361   CHKERRQ(DMGetLocalVector(da,&xlocal));
362   /*
363      Scatter ghost points to local vector, using the 2-step process
364         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
365      By placing code between these two statements, computations can
366      be done while messages are in transition.
367   */
368   CHKERRQ(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal));
369   CHKERRQ(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal));
370 
371   /*
372      Get pointers to vector data.
373        - The vector xlocal includes ghost point; the vectors x and f do
374          NOT include ghost points.
375        - Using DMDAVecGetArray() allows accessing the values using global ordering
376   */
377   CHKERRQ(DMDAVecGetArray(da,xlocal,&xx));
378   CHKERRQ(DMDAVecGetArray(da,f,&ff));
379   CHKERRQ(DMDAVecGetArray(da,user->F,&FF));
380 
381   /*
382      Get local grid boundaries (for 1-dimensional DMDA):
383        xs, xm  - starting grid index, width of local grid (no ghost points)
384   */
385   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
386   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
387 
388   /*
389      Set function values for boundary points; define local interior grid point range:
390         xsi - starting interior grid index
391         xei - ending interior grid index
392   */
393   if (xs == 0) { /* left boundary */
394     ff[0] = xx[0];
395     xs++;xm--;
396   }
397   if (xs+xm == M) {  /* right boundary */
398     ff[xs+xm-1] = xx[xs+xm-1] - 1.0;
399     xm--;
400   }
401 
402   /*
403      Compute function over locally owned part of the grid (interior points only)
404   */
405   d = 1.0/(user->h*user->h);
406   for (i=xs; i<xs+xm; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
407 
408   /*
409      Restore vectors
410   */
411   CHKERRQ(DMDAVecRestoreArray(da,xlocal,&xx));
412   CHKERRQ(DMDAVecRestoreArray(da,f,&ff));
413   CHKERRQ(DMDAVecRestoreArray(da,user->F,&FF));
414   CHKERRQ(DMRestoreLocalVector(da,&xlocal));
415   PetscFunctionReturn(0);
416 }
417 
418 /* ------------------------------------------------------------------- */
419 /*
420    FormJacobian - Evaluates Jacobian matrix.
421 
422    Input Parameters:
423 .  snes - the SNES context
424 .  x - input vector
425 .  dummy - optional user-defined context (not used here)
426 
427    Output Parameters:
428 .  jac - Jacobian matrix
429 .  B - optionally different preconditioning matrix
430 .  flag - flag indicating matrix structure
431 */
432 PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
433 {
434   ApplicationCtx *user = (ApplicationCtx*) ctx;
435   PetscScalar    *xx,d,A[3];
436   PetscInt       i,j[3],M,xs,xm;
437   DM             da = user->da;
438 
439   PetscFunctionBeginUser;
440   if (user->sjerr) {
441     CHKERRQ(SNESSetJacobianDomainError(snes));
442     PetscFunctionReturn(0);
443   }
444   /*
445      Get pointer to vector data
446   */
447   CHKERRQ(DMDAVecGetArrayRead(da,x,&xx));
448   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
449 
450   /*
451     Get range of locally owned matrix
452   */
453   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
454 
455   /*
456      Determine starting and ending local indices for interior grid points.
457      Set Jacobian entries for boundary points.
458   */
459 
460   if (xs == 0) {  /* left boundary */
461     i = 0; A[0] = 1.0;
462 
463     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
464     xs++;xm--;
465   }
466   if (xs+xm == M) { /* right boundary */
467     i    = M-1;
468     A[0] = 1.0;
469     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
470     xm--;
471   }
472 
473   /*
474      Interior grid points
475       - Note that in this case we set all elements for a particular
476         row at once.
477   */
478   d = 1.0/(user->h*user->h);
479   for (i=xs; i<xs+xm; i++) {
480     j[0] = i - 1; j[1] = i; j[2] = i + 1;
481     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
482     CHKERRQ(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
483   }
484 
485   /*
486      Assemble matrix, using the 2-step process:
487        MatAssemblyBegin(), MatAssemblyEnd().
488      By placing code between these two statements, computations can be
489      done while messages are in transition.
490 
491      Also, restore vector.
492   */
493 
494   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
495   CHKERRQ(DMDAVecRestoreArrayRead(da,x,&xx));
496   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
497 
498   PetscFunctionReturn(0);
499 }
500 
501 /* ------------------------------------------------------------------- */
502 /*
503    Monitor - Optional user-defined monitoring routine that views the
504    current iterate with an x-window plot. Set by SNESMonitorSet().
505 
506    Input Parameters:
507    snes - the SNES context
508    its - iteration number
509    norm - 2-norm function value (may be estimated)
510    ctx - optional user-defined context for private data for the
511          monitor routine, as set by SNESMonitorSet()
512 
513    Note:
514    See the manpage for PetscViewerDrawOpen() for useful runtime options,
515    such as -nox to deactivate all x-window output.
516  */
517 PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
518 {
519   MonitorCtx     *monP = (MonitorCtx*) ctx;
520   Vec            x;
521 
522   PetscFunctionBeginUser;
523   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,(double)fnorm));
524   CHKERRQ(SNESGetSolution(snes,&x));
525   CHKERRQ(VecView(x,monP->viewer));
526   PetscFunctionReturn(0);
527 }
528 
529 /* ------------------------------------------------------------------- */
530 /*
531    PreCheck - Optional user-defined routine that checks the validity of
532    candidate steps of a line search method.  Set by SNESLineSearchSetPreCheck().
533 
534    Input Parameters:
535    snes - the SNES context
536    xcurrent - current solution
537    y - search direction and length
538 
539    Output Parameters:
540    y         - proposed step (search direction and length) (possibly changed)
541    changed_y - tells if the step has changed or not
542  */
543 PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx)
544 {
545   PetscFunctionBeginUser;
546   *changed_y = PETSC_FALSE;
547   PetscFunctionReturn(0);
548 }
549 
550 /* ------------------------------------------------------------------- */
551 /*
552    PostCheck - Optional user-defined routine that checks the validity of
553    candidate steps of a line search method.  Set by SNESLineSearchSetPostCheck().
554 
555    Input Parameters:
556    snes - the SNES context
557    ctx  - optional user-defined context for private data for the
558           monitor routine, as set by SNESLineSearchSetPostCheck()
559    xcurrent - current solution
560    y - search direction and length
561    x    - the new candidate iterate
562 
563    Output Parameters:
564    y    - proposed step (search direction and length) (possibly changed)
565    x    - current iterate (possibly modified)
566 
567  */
568 PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
569 {
570   PetscInt       i,iter,xs,xm;
571   StepCheckCtx   *check;
572   ApplicationCtx *user;
573   PetscScalar    *xa,*xa_last,tmp;
574   PetscReal      rdiff;
575   DM             da;
576   SNES           snes;
577 
578   PetscFunctionBeginUser;
579   *changed_x = PETSC_FALSE;
580   *changed_y = PETSC_FALSE;
581 
582   CHKERRQ(SNESLineSearchGetSNES(linesearch, &snes));
583   check = (StepCheckCtx*)ctx;
584   user  = check->user;
585   CHKERRQ(SNESGetIterationNumber(snes,&iter));
586 
587   /* iteration 1 indicates we are working on the second iteration */
588   if (iter > 0) {
589     da   = user->da;
590     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",iter,(double)check->tolerance));
591 
592     /* Access local array data */
593     CHKERRQ(DMDAVecGetArray(da,check->last_step,&xa_last));
594     CHKERRQ(DMDAVecGetArray(da,x,&xa));
595     CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
596 
597     /*
598        If we fail the user-defined check for validity of the candidate iterate,
599        then modify the iterate as we like.  (Note that the particular modification
600        below is intended simply to demonstrate how to manipulate this data, not
601        as a meaningful or appropriate choice.)
602     */
603     for (i=xs; i<xs+xm; i++) {
604       if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance;
605       else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]);
606       if (rdiff > check->tolerance) {
607         tmp        = xa[i];
608         xa[i]      = .5*(xa[i] + xa_last[i]);
609         *changed_x = PETSC_TRUE;
610         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"  Altering entry %D: x=%g, x_last=%g, diff=%g, x_new=%g\n",i,(double)PetscAbsScalar(tmp),(double)PetscAbsScalar(xa_last[i]),(double)rdiff,(double)PetscAbsScalar(xa[i])));
611       }
612     }
613     CHKERRQ(DMDAVecRestoreArray(da,check->last_step,&xa_last));
614     CHKERRQ(DMDAVecRestoreArray(da,x,&xa));
615   }
616   CHKERRQ(VecCopy(x,check->last_step));
617   PetscFunctionReturn(0);
618 }
619 
620 /* ------------------------------------------------------------------- */
621 /*
622    PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used
623    e.g,
624      mpiexec -n 8 ./ex3 -nox -n 10000 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 4 -sub_ksp_type gmres -sub_ksp_max_it 3 -post_setsubksp -sub_ksp_rtol 1.e-16
625    Set by SNESLineSearchSetPostCheck().
626 
627    Input Parameters:
628    linesearch - the LineSearch context
629    xcurrent - current solution
630    y - search direction and length
631    x    - the new candidate iterate
632 
633    Output Parameters:
634    y    - proposed step (search direction and length) (possibly changed)
635    x    - current iterate (possibly modified)
636 
637  */
638 PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
639 {
640   SetSubKSPCtx   *check;
641   PetscInt       iter,its,sub_its,maxit;
642   KSP            ksp,sub_ksp,*sub_ksps;
643   PC             pc;
644   PetscReal      ksp_ratio;
645   SNES           snes;
646 
647   PetscFunctionBeginUser;
648   CHKERRQ(SNESLineSearchGetSNES(linesearch, &snes));
649   check   = (SetSubKSPCtx*)ctx;
650   CHKERRQ(SNESGetIterationNumber(snes,&iter));
651   CHKERRQ(SNESGetKSP(snes,&ksp));
652   CHKERRQ(KSPGetPC(ksp,&pc));
653   CHKERRQ(PCBJacobiGetSubKSP(pc,NULL,NULL,&sub_ksps));
654   sub_ksp = sub_ksps[0];
655   CHKERRQ(KSPGetIterationNumber(ksp,&its));      /* outer KSP iteration number */
656   CHKERRQ(KSPGetIterationNumber(sub_ksp,&sub_its)); /* inner KSP iteration number */
657 
658   if (iter) {
659     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"    ...PostCheck snes iteration %D, ksp_it %D %D, subksp_it %D\n",iter,check->its0,its,sub_its));
660     ksp_ratio = ((PetscReal)(its))/check->its0;
661     maxit     = (PetscInt)(ksp_ratio*sub_its + 0.5);
662     if (maxit < 2) maxit = 2;
663     CHKERRQ(KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit));
664     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"    ...ksp_ratio %g, new maxit %D\n\n",(double)ksp_ratio,maxit));
665   }
666   check->its0 = its; /* save current outer KSP iteration number */
667   PetscFunctionReturn(0);
668 }
669 
670 /* ------------------------------------------------------------------- */
671 /*
672    MatrixFreePreconditioner - This routine demonstrates the use of a
673    user-provided preconditioner.  This code implements just the null
674    preconditioner, which of course is not recommended for general use.
675 
676    Input Parameters:
677 +  pc - preconditioner
678 -  x - input vector
679 
680    Output Parameter:
681 .  y - preconditioned vector
682 */
683 PetscErrorCode MatrixFreePreconditioner(PC pc,Vec x,Vec y)
684 {
685   CHKERRQ(VecCopy(x,y));
686   return 0;
687 }
688 
689 /*TEST
690 
691    test:
692       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
693 
694    test:
695       suffix: 2
696       nsize: 3
697       args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
698 
699    test:
700       suffix: 3
701       nsize: 2
702       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
703 
704    test:
705       suffix: 4
706       args: -nox -pre_check_iterates -post_check_iterates
707 
708    test:
709       suffix: 5
710       requires: double !complex !single
711       nsize: 2
712       args: -nox -snes_test_jacobian  -snes_test_jacobian_view
713 
714    test:
715       suffix: 6
716       requires: double !complex !single
717       nsize: 4
718       args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1
719 
720    test:
721       suffix: 7
722       requires: double !complex !single
723       nsize: 4
724       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1
725 
726    test:
727       suffix: 8
728       requires: double !complex !single
729       nsize: 4
730       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1
731 
732    test:
733       suffix: 9
734       requires: double !complex !single
735       nsize: 4
736       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1
737 
738    test:
739       suffix: 10
740       requires: double !complex !single
741       nsize: 4
742       args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1
743 
744    test:
745       suffix: 11
746       requires: double !complex !single
747       nsize: 4
748       args: -test_jacobian_domain_error -snes_converged_reason -snes_type ms -snes_ms_type m62 -snes_ms_damping 0.9 -snes_check_jacobian_domain_error 1
749 
750    test:
751       suffix: 12
752       args: -view_initial
753       filter: grep -v "type:"
754 
755    test:
756       suffix: 13
757       requires: double !complex !single
758       nsize: 4
759       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1
760 
761 TEST*/
762