Actual source code: ex3.c
petsc-3.8.0 2017-09-26
2: static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel.\n\
3: This example employs a user-defined monitoring routine and optionally a user-defined\n\
4: routine to check candidate iterates produced by line search routines. This code also\n\
5: The command line options include:\n\
6: -pre_check_iterates : activate checking of iterates\n\
7: -post_check_iterates : activate checking of iterates\n\
8: -check_tol <tol>: set tolerance for iterate checking\n\n";
10: /*T
11: Concepts: SNES^basic parallel example
12: Concepts: SNES^setting a user-defined monitoring routine
13: Processors: n
14: T*/
16: /*
17: Include "petscdraw.h" so that we can use distributed arrays (DMDAs).
18: Include "petscdraw.h" so that we can use PETSc drawing routines.
19: Include "petscsnes.h" so that we can use SNES solvers. Note that this
20: file automatically includes:
21: petscsys.h - base PETSc routines petscvec.h - vectors
22: petscmat.h - matrices
23: petscis.h - index sets petscksp.h - Krylov subspace methods
24: petscviewer.h - viewers petscpc.h - preconditioners
25: petscksp.h - linear solvers
26: */
28: #include <petscdm.h>
29: #include <petscdmda.h>
30: #include <petscsnes.h>
32: /*
33: User-defined routines. Note that immediately before each routine below,
34: If defined, this macro is used in the PETSc error handlers to provide a
35: complete traceback of routine names. All PETSc library routines use this
36: macro, and users can optionally employ it as well in their application
37: codes. Note that users can get a traceback of PETSc errors regardless of
38: provides the added traceback detail of the application routine names.
39: */
40: PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
41: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
42: PetscErrorCode FormInitialGuess(Vec);
43: PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
44: PetscErrorCode PreCheck(SNESLineSearch,Vec,Vec,PetscBool*,void*);
45: PetscErrorCode PostCheck(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
46: PetscErrorCode PostSetSubKSP(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
47: PetscErrorCode MatrixFreePreconditioner(PC,Vec,Vec);
49: /*
50: User-defined application context
51: */
52: typedef struct {
53: DM da; /* distributed array */
54: Vec F; /* right-hand-side of PDE */
55: PetscMPIInt rank; /* rank of processor */
56: PetscMPIInt size; /* size of communicator */
57: PetscReal h; /* mesh spacing */
58: } ApplicationCtx;
60: /*
61: User-defined context for monitoring
62: */
63: typedef struct {
64: PetscViewer viewer;
65: } MonitorCtx;
67: /*
68: User-defined context for checking candidate iterates that are
69: determined by line search methods
70: */
71: typedef struct {
72: Vec last_step; /* previous iterate */
73: PetscReal tolerance; /* tolerance for changes between successive iterates */
74: ApplicationCtx *user;
75: } StepCheckCtx;
77: typedef struct {
78: PetscInt its0; /* num of prevous outer KSP iterations */
79: } SetSubKSPCtx;
81: int main(int argc,char **argv)
82: {
83: SNES snes; /* SNES context */
84: SNESLineSearch linesearch; /* SNESLineSearch context */
85: Mat J; /* Jacobian matrix */
86: ApplicationCtx ctx; /* user-defined context */
87: Vec x,r,U,F; /* vectors */
88: MonitorCtx monP; /* monitoring context */
89: StepCheckCtx checkP; /* step-checking context */
90: SetSubKSPCtx checkP1;
91: PetscBool pre_check,post_check,post_setsubksp; /* flag indicating whether we're checking candidate iterates */
92: PetscScalar xp,*FF,*UU,none = -1.0;
94: PetscInt its,N = 5,i,maxit,maxf,xs,xm;
95: PetscReal abstol,rtol,stol,norm;
96: PetscBool flg;
99: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
100: MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank);
101: MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size);
102: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
103: ctx.h = 1.0/(N-1);
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Create nonlinear solver context
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: SNESCreate(PETSC_COMM_WORLD,&snes);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Create vector data structures; set function evaluation routine
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: /*
116: Create distributed array (DMDA) to manage parallel grid and vectors
117: */
118: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da);
119: DMSetFromOptions(ctx.da);
120: DMSetUp(ctx.da);
122: /*
123: Extract global and local vectors from DMDA; then duplicate for remaining
124: vectors that are the same types
125: */
126: DMCreateGlobalVector(ctx.da,&x);
127: VecDuplicate(x,&r);
128: VecDuplicate(x,&F); ctx.F = F;
129: VecDuplicate(x,&U);
131: /*
132: Set function evaluation routine and vector. Whenever the nonlinear
133: solver needs to compute the nonlinear function, it will call this
134: routine.
135: - Note that the final routine argument is the user-defined
136: context that provides application-specific data for the
137: function evaluation routine.
138: */
139: SNESSetFunction(snes,r,FormFunction,&ctx);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Create matrix data structure; set Jacobian evaluation routine
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: MatCreate(PETSC_COMM_WORLD,&J);
146: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N);
147: MatSetFromOptions(J);
148: MatSeqAIJSetPreallocation(J,3,NULL);
149: MatMPIAIJSetPreallocation(J,3,NULL,3,NULL);
151: /*
152: Set Jacobian matrix data structure and default Jacobian evaluation
153: routine. Whenever the nonlinear solver needs to compute the
154: Jacobian matrix, it will call this routine.
155: - Note that the final routine argument is the user-defined
156: context that provides application-specific data for the
157: Jacobian evaluation routine.
158: */
159: SNESSetJacobian(snes,J,J,FormJacobian,&ctx);
161: /*
162: Optional allow user provided preconditioner
163: */
164: PetscOptionsHasName(NULL,NULL,"-user_precond",&flg);
165: if (flg) {
166: KSP ksp;
167: PC pc;
168: SNESGetKSP(snes,&ksp);
169: KSPGetPC(ksp,&pc);
170: PCSetType(pc,PCSHELL);
171: PCShellSetApply(pc,MatrixFreePreconditioner);
172: }
174: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: Customize nonlinear solver; set runtime options
176: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178: /*
179: Set an optional user-defined monitoring routine
180: */
181: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);
182: SNESMonitorSet(snes,Monitor,&monP,0);
184: /*
185: Set names for some vectors to facilitate monitoring (optional)
186: */
187: PetscObjectSetName((PetscObject)x,"Approximate Solution");
188: PetscObjectSetName((PetscObject)U,"Exact Solution");
190: /*
191: Set SNES/KSP/KSP/PC runtime options, e.g.,
192: -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
193: */
194: SNESSetFromOptions(snes);
196: /*
197: Set an optional user-defined routine to check the validity of candidate
198: iterates that are determined by line search methods
199: */
200: SNESGetLineSearch(snes, &linesearch);
201: PetscOptionsHasName(NULL,NULL,"-post_check_iterates",&post_check);
203: if (post_check) {
204: PetscPrintf(PETSC_COMM_WORLD,"Activating post step checking routine\n");
205: SNESLineSearchSetPostCheck(linesearch,PostCheck,&checkP);
206: VecDuplicate(x,&(checkP.last_step));
208: checkP.tolerance = 1.0;
209: checkP.user = &ctx;
211: PetscOptionsGetReal(NULL,NULL,"-check_tol",&checkP.tolerance,NULL);
212: }
214: PetscOptionsHasName(NULL,NULL,"-post_setsubksp",&post_setsubksp);
215: if (post_setsubksp) {
216: PetscPrintf(PETSC_COMM_WORLD,"Activating post setsubksp\n");
217: SNESLineSearchSetPostCheck(linesearch,PostSetSubKSP,&checkP1);
218: }
220: PetscOptionsHasName(NULL,NULL,"-pre_check_iterates",&pre_check);
221: if (pre_check) {
222: PetscPrintf(PETSC_COMM_WORLD,"Activating pre step checking routine\n");
223: SNESLineSearchSetPreCheck(linesearch,PreCheck,&checkP);
224: }
226: /*
227: Print parameters used for convergence testing (optional) ... just
228: to demonstrate this routine; this information is also printed with
229: the option -snes_view
230: */
231: SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);
232: 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: Initialize application:
236: Store right-hand-side of PDE and exact solution
237: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239: /*
240: Get local grid boundaries (for 1-dimensional DMDA):
241: xs, xm - starting grid index, width of local grid (no ghost points)
242: */
243: DMDAGetCorners(ctx.da,&xs,NULL,NULL,&xm,NULL,NULL);
245: /*
246: Get pointers to vector data
247: */
248: DMDAVecGetArray(ctx.da,F,&FF);
249: DMDAVecGetArray(ctx.da,U,&UU);
251: /*
252: Compute local vector entries
253: */
254: xp = ctx.h*xs;
255: for (i=xs; i<xs+xm; i++) {
256: FF[i] = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
257: UU[i] = xp*xp*xp;
258: xp += ctx.h;
259: }
261: /*
262: Restore vectors
263: */
264: DMDAVecRestoreArray(ctx.da,F,&FF);
265: DMDAVecRestoreArray(ctx.da,U,&UU);
267: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268: Evaluate initial guess; then solve nonlinear system
269: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271: /*
272: Note: The user should initialize the vector, x, with the initial guess
273: for the nonlinear solver prior to calling SNESSolve(). In particular,
274: to employ an initial guess of zero, the user should explicitly set
275: this vector to zero by calling VecSet().
276: */
277: FormInitialGuess(x);
278: SNESSolve(snes,NULL,x);
279: SNESGetIterationNumber(snes,&its);
280: PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);
282: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283: Check solution and clean up
284: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
286: /*
287: Check the error
288: */
289: VecAXPY(x,none,U);
290: VecNorm(x,NORM_2,&norm);
291: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its);
293: /*
294: Free work space. All PETSc objects should be destroyed when they
295: are no longer needed.
296: */
297: PetscViewerDestroy(&monP.viewer);
298: if (post_check) {VecDestroy(&checkP.last_step);}
299: VecDestroy(&x);
300: VecDestroy(&r);
301: VecDestroy(&U);
302: VecDestroy(&F);
303: MatDestroy(&J);
304: SNESDestroy(&snes);
305: DMDestroy(&ctx.da);
306: PetscFinalize();
307: return ierr;
308: }
309: /* ------------------------------------------------------------------- */
310: /*
311: FormInitialGuess - Computes initial guess.
313: Input/Output Parameter:
314: . x - the solution vector
315: */
316: PetscErrorCode FormInitialGuess(Vec x)
317: {
319: PetscScalar pfive = .50;
322: VecSet(x,pfive);
323: return(0);
324: }
325: /* ------------------------------------------------------------------- */
326: /*
327: FormFunction - Evaluates nonlinear function, F(x).
329: Input Parameters:
330: . snes - the SNES context
331: . x - input vector
332: . ctx - optional user-defined context, as set by SNESSetFunction()
334: Output Parameter:
335: . f - function vector
337: Note:
338: The user-defined context can contain any application-specific
339: data needed for the function evaluation.
340: */
341: PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
342: {
343: ApplicationCtx *user = (ApplicationCtx*) ctx;
344: DM da = user->da;
345: PetscScalar *xx,*ff,*FF,d;
347: PetscInt i,M,xs,xm;
348: Vec xlocal;
351: DMGetLocalVector(da,&xlocal);
352: /*
353: Scatter ghost points to local vector, using the 2-step process
354: DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
355: By placing code between these two statements, computations can
356: be done while messages are in transition.
357: */
358: DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);
359: DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);
361: /*
362: Get pointers to vector data.
363: - The vector xlocal includes ghost point; the vectors x and f do
364: NOT include ghost points.
365: - Using DMDAVecGetArray() allows accessing the values using global ordering
366: */
367: DMDAVecGetArray(da,xlocal,&xx);
368: DMDAVecGetArray(da,f,&ff);
369: DMDAVecGetArray(da,user->F,&FF);
371: /*
372: Get local grid boundaries (for 1-dimensional DMDA):
373: xs, xm - starting grid index, width of local grid (no ghost points)
374: */
375: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
376: DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
378: /*
379: Set function values for boundary points; define local interior grid point range:
380: xsi - starting interior grid index
381: xei - ending interior grid index
382: */
383: if (xs == 0) { /* left boundary */
384: ff[0] = xx[0];
385: xs++;xm--;
386: }
387: if (xs+xm == M) { /* right boundary */
388: ff[xs+xm-1] = xx[xs+xm-1] - 1.0;
389: xm--;
390: }
392: /*
393: Compute function over locally owned part of the grid (interior points only)
394: */
395: d = 1.0/(user->h*user->h);
396: 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];
398: /*
399: Restore vectors
400: */
401: DMDAVecRestoreArray(da,xlocal,&xx);
402: DMDAVecRestoreArray(da,f,&ff);
403: DMDAVecRestoreArray(da,user->F,&FF);
404: DMRestoreLocalVector(da,&xlocal);
405: return(0);
406: }
407: /* ------------------------------------------------------------------- */
408: /*
409: FormJacobian - Evaluates Jacobian matrix.
411: Input Parameters:
412: . snes - the SNES context
413: . x - input vector
414: . dummy - optional user-defined context (not used here)
416: Output Parameters:
417: . jac - Jacobian matrix
418: . B - optionally different preconditioning matrix
419: . flag - flag indicating matrix structure
420: */
421: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
422: {
423: ApplicationCtx *user = (ApplicationCtx*) ctx;
424: PetscScalar *xx,d,A[3];
426: PetscInt i,j[3],M,xs,xm;
427: DM da = user->da;
430: /*
431: Get pointer to vector data
432: */
433: DMDAVecGetArrayRead(da,x,&xx);
434: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
436: /*
437: Get range of locally owned matrix
438: */
439: DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
441: /*
442: Determine starting and ending local indices for interior grid points.
443: Set Jacobian entries for boundary points.
444: */
446: if (xs == 0) { /* left boundary */
447: i = 0; A[0] = 1.0;
449: MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES);
450: xs++;xm--;
451: }
452: if (xs+xm == M) { /* right boundary */
453: i = M-1;
454: A[0] = 1.0;
455: MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES);
456: xm--;
457: }
459: /*
460: Interior grid points
461: - Note that in this case we set all elements for a particular
462: row at once.
463: */
464: d = 1.0/(user->h*user->h);
465: for (i=xs; i<xs+xm; i++) {
466: j[0] = i - 1; j[1] = i; j[2] = i + 1;
467: A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
468: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
469: }
471: /*
472: Assemble matrix, using the 2-step process:
473: MatAssemblyBegin(), MatAssemblyEnd().
474: By placing code between these two statements, computations can be
475: done while messages are in transition.
477: Also, restore vector.
478: */
480: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
481: DMDAVecRestoreArrayRead(da,x,&xx);
482: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
484: return(0);
485: }
486: /* ------------------------------------------------------------------- */
487: /*
488: Monitor - Optional user-defined monitoring routine that views the
489: current iterate with an x-window plot. Set by SNESMonitorSet().
491: Input Parameters:
492: snes - the SNES context
493: its - iteration number
494: norm - 2-norm function value (may be estimated)
495: ctx - optional user-defined context for private data for the
496: monitor routine, as set by SNESMonitorSet()
498: Note:
499: See the manpage for PetscViewerDrawOpen() for useful runtime options,
500: such as -nox to deactivate all x-window output.
501: */
502: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
503: {
505: MonitorCtx *monP = (MonitorCtx*) ctx;
506: Vec x;
509: PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,(double)fnorm);
510: SNESGetSolution(snes,&x);
511: VecView(x,monP->viewer);
512: return(0);
513: }
515: /* ------------------------------------------------------------------- */
516: /*
517: PreCheck - Optional user-defined routine that checks the validity of
518: candidate steps of a line search method. Set by SNESLineSearchSetPreCheck().
520: Input Parameters:
521: snes - the SNES context
522: xcurrent - current solution
523: y - search direction and length
525: Output Parameters:
526: y - proposed step (search direction and length) (possibly changed)
527: changed_y - tells if the step has changed or not
528: */
529: PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx)
530: {
532: *changed_y = PETSC_FALSE;
533: return(0);
534: }
536: /* ------------------------------------------------------------------- */
537: /*
538: PostCheck - Optional user-defined routine that checks the validity of
539: candidate steps of a line search method. Set by SNESLineSearchSetPostCheck().
541: Input Parameters:
542: snes - the SNES context
543: ctx - optional user-defined context for private data for the
544: monitor routine, as set by SNESLineSearchSetPostCheck()
545: xcurrent - current solution
546: y - search direction and length
547: x - the new candidate iterate
549: Output Parameters:
550: y - proposed step (search direction and length) (possibly changed)
551: x - current iterate (possibly modified)
553: */
554: PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool *changed_y,PetscBool *changed_x, void * ctx)
555: {
557: PetscInt i,iter,xs,xm;
558: StepCheckCtx *check;
559: ApplicationCtx *user;
560: PetscScalar *xa,*xa_last,tmp;
561: PetscReal rdiff;
562: DM da;
563: SNES snes;
566: *changed_x = PETSC_FALSE;
567: *changed_y = PETSC_FALSE;
569: SNESLineSearchGetSNES(linesearch, &snes);
570: check = (StepCheckCtx*)ctx;
571: user = check->user;
572: SNESGetIterationNumber(snes,&iter);
574: /* iteration 1 indicates we are working on the second iteration */
575: if (iter > 0) {
576: da = user->da;
577: PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",iter,(double)check->tolerance);
579: /* Access local array data */
580: DMDAVecGetArray(da,check->last_step,&xa_last);
581: DMDAVecGetArray(da,x,&xa);
582: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
584: /*
585: If we fail the user-defined check for validity of the candidate iterate,
586: then modify the iterate as we like. (Note that the particular modification
587: below is intended simply to demonstrate how to manipulate this data, not
588: as a meaningful or appropriate choice.)
589: */
590: for (i=xs; i<xs+xm; i++) {
591: if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance;
592: else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]);
593: if (rdiff > check->tolerance) {
594: tmp = xa[i];
595: xa[i] = .5*(xa[i] + xa_last[i]);
596: *changed_x = PETSC_TRUE;
597: 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]));
598: }
599: }
600: DMDAVecRestoreArray(da,check->last_step,&xa_last);
601: DMDAVecRestoreArray(da,x,&xa);
602: }
603: VecCopy(x,check->last_step);
604: return(0);
605: }
608: /* ------------------------------------------------------------------- */
609: /*
610: PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used
611: e.g,
612: 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
613: Set by SNESLineSearchSetPostCheck().
615: Input Parameters:
616: linesearch - the LineSearch context
617: xcurrent - current solution
618: y - search direction and length
619: x - the new candidate iterate
621: Output Parameters:
622: y - proposed step (search direction and length) (possibly changed)
623: x - current iterate (possibly modified)
625: */
626: PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool *changed_y,PetscBool *changed_x, void * ctx)
627: {
629: SetSubKSPCtx *check;
630: PetscInt iter,its,sub_its,maxit;
631: KSP ksp,sub_ksp,*sub_ksps;
632: PC pc;
633: PetscReal ksp_ratio;
634: SNES snes;
637: SNESLineSearchGetSNES(linesearch, &snes);
638: check = (SetSubKSPCtx*)ctx;
639: SNESGetIterationNumber(snes,&iter);
640: SNESGetKSP(snes,&ksp);
641: KSPGetPC(ksp,&pc);
642: PCBJacobiGetSubKSP(pc,NULL,NULL,&sub_ksps);
643: sub_ksp = sub_ksps[0];
644: KSPGetIterationNumber(ksp,&its); /* outer KSP iteration number */
645: KSPGetIterationNumber(sub_ksp,&sub_its); /* inner KSP iteration number */
647: if (iter) {
648: PetscPrintf(PETSC_COMM_WORLD," ...PostCheck snes iteration %D, ksp_it %d %d, subksp_it %d\n",iter,check->its0,its,sub_its);
649: ksp_ratio = ((PetscReal)(its))/check->its0;
650: maxit = (PetscInt)(ksp_ratio*sub_its + 0.5);
651: if (maxit < 2) maxit = 2;
652: KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit);
653: PetscPrintf(PETSC_COMM_WORLD," ...ksp_ratio %g, new maxit %d\n\n",ksp_ratio,maxit);
654: }
655: check->its0 = its; /* save current outer KSP iteration number */
656: return(0);
657: }
659: /* ------------------------------------------------------------------- */
660: /*
661: MatrixFreePreconditioner - This routine demonstrates the use of a
662: user-provided preconditioner. This code implements just the null
663: preconditioner, which of course is not recommended for general use.
665: Input Parameters:
666: + pc - preconditioner
667: - x - input vector
669: Output Parameter:
670: . y - preconditioned vector
671: */
672: PetscErrorCode MatrixFreePreconditioner(PC pc,Vec x,Vec y)
673: {
675: VecCopy(x,y);
676: return 0;
677: }