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