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