xref: /petsc/src/snes/tutorials/ex3.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
98*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank));
99*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size));
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
101c4762a1bSJed Brown   ctx.h = 1.0/(N-1);
102c4762a1bSJed Brown   ctx.sjerr = PETSC_FALSE;
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_jacobian_domain_error",&ctx.sjerr,NULL));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Create nonlinear solver context
108c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109c4762a1bSJed Brown 
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
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   */
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(ctx.da));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(ctx.da));
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   */
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(ctx.da,&x));
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F)); ctx.F = F;
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
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   */
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormFunction,&ctx));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
144c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145c4762a1bSJed Brown 
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(J,3,NULL,3,NULL));
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   */
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,&ctx));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /*
163c4762a1bSJed Brown      Optionally allow user-provided preconditioner
164c4762a1bSJed Brown    */
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-user_precond",&flg));
166c4762a1bSJed Brown   if (flg) {
167c4762a1bSJed Brown     KSP ksp;
168c4762a1bSJed Brown     PC  pc;
169*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetKSP(snes,&ksp));
170*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
171*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCSHELL));
172*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellSetApply(pc,MatrixFreePreconditioner));
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   */
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorSet(snes,Monitor,&monP,0));
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /*
186c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
187c4762a1bSJed Brown   */
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
189*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution"));
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   */
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
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   */
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLineSearch(snes, &linesearch));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-post_check_iterates",&post_check));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   if (post_check) {
205*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating post step checking routine\n"));
206*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESLineSearchSetPostCheck(linesearch,PostCheck,&checkP));
207*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&(checkP.last_step)));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown     checkP.tolerance = 1.0;
210c4762a1bSJed Brown     checkP.user      = &ctx;
211c4762a1bSJed Brown 
212*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-check_tol",&checkP.tolerance,NULL));
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown 
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-post_setsubksp",&post_setsubksp));
216c4762a1bSJed Brown   if (post_setsubksp) {
217*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating post setsubksp\n"));
218*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESLineSearchSetPostCheck(linesearch,PostSetSubKSP,&checkP1));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-pre_check_iterates",&pre_check));
222c4762a1bSJed Brown   if (pre_check) {
223*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating pre step checking routine\n"));
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESLineSearchSetPreCheck(linesearch,PreCheck,&checkP));
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   */
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
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   */
244*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(ctx.da,&xs,NULL,NULL,&xm,NULL,NULL));
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /*
247c4762a1bSJed Brown      Get pointers to vector data
248c4762a1bSJed Brown   */
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(ctx.da,F,&FF));
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(ctx.da,U,&UU));
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   */
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(ctx.da,F,&FF));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(ctx.da,U,&UU));
267c0558f20SBarry Smith   if (viewinitial) {
268*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(U,0));
269*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(F,0));
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   */
282*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288c4762a1bSJed Brown      Check solution and clean up
289c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   /*
292c4762a1bSJed Brown      Check the error
293c4762a1bSJed Brown   */
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,none,U));
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
296*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its));
297c4762a1bSJed Brown   if (ctx.sjerr) {
298c4762a1bSJed Brown     SNESType       snestype;
299*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetType(snes,&snestype));
300*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"SNES Type %s\n",snestype));
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   */
307*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&monP.viewer));
308*5f80ce2aSJacob Faibussowitsch   if (post_check) CHKERRQ(VecDestroy(&checkP.last_step));
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&F));
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
315*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&ctx.da));
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   PetscScalar    pfive = .50;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   PetscFunctionBeginUser;
332*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
333c4762a1bSJed Brown   PetscFunctionReturn(0);
334c4762a1bSJed Brown }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown /* ------------------------------------------------------------------- */
337c4762a1bSJed Brown /*
338c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    Input Parameters:
341c4762a1bSJed Brown .  snes - the SNES context
342c4762a1bSJed Brown .  x - input vector
343c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
344c4762a1bSJed Brown 
345c4762a1bSJed Brown    Output Parameter:
346c4762a1bSJed Brown .  f - function vector
347c4762a1bSJed Brown 
348c4762a1bSJed Brown    Note:
349c4762a1bSJed Brown    The user-defined context can contain any application-specific
350c4762a1bSJed Brown    data needed for the function evaluation.
351c4762a1bSJed Brown */
352c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
353c4762a1bSJed Brown {
354c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx*) ctx;
355c4762a1bSJed Brown   DM             da    = user->da;
356c4762a1bSJed Brown   PetscScalar    *xx,*ff,*FF,d;
357c4762a1bSJed Brown   PetscInt       i,M,xs,xm;
358c4762a1bSJed Brown   Vec            xlocal;
359c4762a1bSJed Brown 
360c4762a1bSJed Brown   PetscFunctionBeginUser;
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&xlocal));
362c4762a1bSJed Brown   /*
363c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
364c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
365c4762a1bSJed Brown      By placing code between these two statements, computations can
366c4762a1bSJed Brown      be done while messages are in transition.
367c4762a1bSJed Brown   */
368*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal));
369*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /*
372c4762a1bSJed Brown      Get pointers to vector data.
373c4762a1bSJed Brown        - The vector xlocal includes ghost point; the vectors x and f do
374c4762a1bSJed Brown          NOT include ghost points.
375c4762a1bSJed Brown        - Using DMDAVecGetArray() allows accessing the values using global ordering
376c4762a1bSJed Brown   */
377*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,xlocal,&xx));
378*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,f,&ff));
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,user->F,&FF));
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   /*
382c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
383c4762a1bSJed Brown        xs, xm  - starting grid index, width of local grid (no ghost points)
384c4762a1bSJed Brown   */
385*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
386*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   /*
389c4762a1bSJed Brown      Set function values for boundary points; define local interior grid point range:
390c4762a1bSJed Brown         xsi - starting interior grid index
391c4762a1bSJed Brown         xei - ending interior grid index
392c4762a1bSJed Brown   */
393c4762a1bSJed Brown   if (xs == 0) { /* left boundary */
394c4762a1bSJed Brown     ff[0] = xx[0];
395c4762a1bSJed Brown     xs++;xm--;
396c4762a1bSJed Brown   }
397c4762a1bSJed Brown   if (xs+xm == M) {  /* right boundary */
398c4762a1bSJed Brown     ff[xs+xm-1] = xx[xs+xm-1] - 1.0;
399c4762a1bSJed Brown     xm--;
400c4762a1bSJed Brown   }
401c4762a1bSJed Brown 
402c4762a1bSJed Brown   /*
403c4762a1bSJed Brown      Compute function over locally owned part of the grid (interior points only)
404c4762a1bSJed Brown   */
405c4762a1bSJed Brown   d = 1.0/(user->h*user->h);
406c4762a1bSJed 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];
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   /*
409c4762a1bSJed Brown      Restore vectors
410c4762a1bSJed Brown   */
411*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,xlocal,&xx));
412*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,f,&ff));
413*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,user->F,&FF));
414*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&xlocal));
415c4762a1bSJed Brown   PetscFunctionReturn(0);
416c4762a1bSJed Brown }
417c4762a1bSJed Brown 
418c4762a1bSJed Brown /* ------------------------------------------------------------------- */
419c4762a1bSJed Brown /*
420c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
421c4762a1bSJed Brown 
422c4762a1bSJed Brown    Input Parameters:
423c4762a1bSJed Brown .  snes - the SNES context
424c4762a1bSJed Brown .  x - input vector
425c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
426c4762a1bSJed Brown 
427c4762a1bSJed Brown    Output Parameters:
428c4762a1bSJed Brown .  jac - Jacobian matrix
429c4762a1bSJed Brown .  B - optionally different preconditioning matrix
430c4762a1bSJed Brown .  flag - flag indicating matrix structure
431c4762a1bSJed Brown */
432c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
433c4762a1bSJed Brown {
434c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx*) ctx;
435c4762a1bSJed Brown   PetscScalar    *xx,d,A[3];
436c4762a1bSJed Brown   PetscInt       i,j[3],M,xs,xm;
437c4762a1bSJed Brown   DM             da = user->da;
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   PetscFunctionBeginUser;
440c4762a1bSJed Brown   if (user->sjerr) {
441*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetJacobianDomainError(snes));
442c4762a1bSJed Brown     PetscFunctionReturn(0);
443c4762a1bSJed Brown   }
444c4762a1bSJed Brown   /*
445c4762a1bSJed Brown      Get pointer to vector data
446c4762a1bSJed Brown   */
447*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,x,&xx));
448*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   /*
451c4762a1bSJed Brown     Get range of locally owned matrix
452c4762a1bSJed Brown   */
453*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   /*
456c4762a1bSJed Brown      Determine starting and ending local indices for interior grid points.
457c4762a1bSJed Brown      Set Jacobian entries for boundary points.
458c4762a1bSJed Brown   */
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   if (xs == 0) {  /* left boundary */
461c4762a1bSJed Brown     i = 0; A[0] = 1.0;
462c4762a1bSJed Brown 
463*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
464c4762a1bSJed Brown     xs++;xm--;
465c4762a1bSJed Brown   }
466c4762a1bSJed Brown   if (xs+xm == M) { /* right boundary */
467c4762a1bSJed Brown     i    = M-1;
468c4762a1bSJed Brown     A[0] = 1.0;
469*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
470c4762a1bSJed Brown     xm--;
471c4762a1bSJed Brown   }
472c4762a1bSJed Brown 
473c4762a1bSJed Brown   /*
474c4762a1bSJed Brown      Interior grid points
475c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
476c4762a1bSJed Brown         row at once.
477c4762a1bSJed Brown   */
478c4762a1bSJed Brown   d = 1.0/(user->h*user->h);
479c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
480c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
481c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
482*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
483c4762a1bSJed Brown   }
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   /*
486c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
487c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
488c4762a1bSJed Brown      By placing code between these two statements, computations can be
489c4762a1bSJed Brown      done while messages are in transition.
490c4762a1bSJed Brown 
491c4762a1bSJed Brown      Also, restore vector.
492c4762a1bSJed Brown   */
493c4762a1bSJed Brown 
494*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
495*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,x,&xx));
496*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   PetscFunctionReturn(0);
499c4762a1bSJed Brown }
500c4762a1bSJed Brown 
501c4762a1bSJed Brown /* ------------------------------------------------------------------- */
502c4762a1bSJed Brown /*
503c4762a1bSJed Brown    Monitor - Optional user-defined monitoring routine that views the
504c4762a1bSJed Brown    current iterate with an x-window plot. Set by SNESMonitorSet().
505c4762a1bSJed Brown 
506c4762a1bSJed Brown    Input Parameters:
507c4762a1bSJed Brown    snes - the SNES context
508c4762a1bSJed Brown    its - iteration number
509c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
510c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
511c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
512c4762a1bSJed Brown 
513c4762a1bSJed Brown    Note:
514c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
515c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
516c4762a1bSJed Brown  */
517c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
518c4762a1bSJed Brown {
519c4762a1bSJed Brown   MonitorCtx     *monP = (MonitorCtx*) ctx;
520c4762a1bSJed Brown   Vec            x;
521c4762a1bSJed Brown 
522c4762a1bSJed Brown   PetscFunctionBeginUser;
523*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,(double)fnorm));
524*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes,&x));
525*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,monP->viewer));
526c4762a1bSJed Brown   PetscFunctionReturn(0);
527c4762a1bSJed Brown }
528c4762a1bSJed Brown 
529c4762a1bSJed Brown /* ------------------------------------------------------------------- */
530c4762a1bSJed Brown /*
531c4762a1bSJed Brown    PreCheck - Optional user-defined routine that checks the validity of
532c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPreCheck().
533c4762a1bSJed Brown 
534c4762a1bSJed Brown    Input Parameters:
535c4762a1bSJed Brown    snes - the SNES context
536c4762a1bSJed Brown    xcurrent - current solution
537c4762a1bSJed Brown    y - search direction and length
538c4762a1bSJed Brown 
539c4762a1bSJed Brown    Output Parameters:
540c4762a1bSJed Brown    y         - proposed step (search direction and length) (possibly changed)
541c4762a1bSJed Brown    changed_y - tells if the step has changed or not
542c4762a1bSJed Brown  */
543c4762a1bSJed Brown PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx)
544c4762a1bSJed Brown {
545c4762a1bSJed Brown   PetscFunctionBeginUser;
546c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
547c4762a1bSJed Brown   PetscFunctionReturn(0);
548c4762a1bSJed Brown }
549c4762a1bSJed Brown 
550c4762a1bSJed Brown /* ------------------------------------------------------------------- */
551c4762a1bSJed Brown /*
552c4762a1bSJed Brown    PostCheck - Optional user-defined routine that checks the validity of
553c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPostCheck().
554c4762a1bSJed Brown 
555c4762a1bSJed Brown    Input Parameters:
556c4762a1bSJed Brown    snes - the SNES context
557c4762a1bSJed Brown    ctx  - optional user-defined context for private data for the
558c4762a1bSJed Brown           monitor routine, as set by SNESLineSearchSetPostCheck()
559c4762a1bSJed Brown    xcurrent - current solution
560c4762a1bSJed Brown    y - search direction and length
561c4762a1bSJed Brown    x    - the new candidate iterate
562c4762a1bSJed Brown 
563c4762a1bSJed Brown    Output Parameters:
564c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
565c4762a1bSJed Brown    x    - current iterate (possibly modified)
566c4762a1bSJed Brown 
567c4762a1bSJed Brown  */
568c4762a1bSJed Brown PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
569c4762a1bSJed Brown {
570c4762a1bSJed Brown   PetscInt       i,iter,xs,xm;
571c4762a1bSJed Brown   StepCheckCtx   *check;
572c4762a1bSJed Brown   ApplicationCtx *user;
573c4762a1bSJed Brown   PetscScalar    *xa,*xa_last,tmp;
574c4762a1bSJed Brown   PetscReal      rdiff;
575c4762a1bSJed Brown   DM             da;
576c4762a1bSJed Brown   SNES           snes;
577c4762a1bSJed Brown 
578c4762a1bSJed Brown   PetscFunctionBeginUser;
579c4762a1bSJed Brown   *changed_x = PETSC_FALSE;
580c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
581c4762a1bSJed Brown 
582*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESLineSearchGetSNES(linesearch, &snes));
583c4762a1bSJed Brown   check = (StepCheckCtx*)ctx;
584c4762a1bSJed Brown   user  = check->user;
585*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&iter));
586c4762a1bSJed Brown 
587c4762a1bSJed Brown   /* iteration 1 indicates we are working on the second iteration */
588c4762a1bSJed Brown   if (iter > 0) {
589c4762a1bSJed Brown     da   = user->da;
590*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",iter,(double)check->tolerance));
591c4762a1bSJed Brown 
592c4762a1bSJed Brown     /* Access local array data */
593*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,check->last_step,&xa_last));
594*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,x,&xa));
595*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
596c4762a1bSJed Brown 
597c4762a1bSJed Brown     /*
598c4762a1bSJed Brown        If we fail the user-defined check for validity of the candidate iterate,
599c4762a1bSJed Brown        then modify the iterate as we like.  (Note that the particular modification
600c4762a1bSJed Brown        below is intended simply to demonstrate how to manipulate this data, not
601c4762a1bSJed Brown        as a meaningful or appropriate choice.)
602c4762a1bSJed Brown     */
603c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
604c4762a1bSJed Brown       if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance;
605c4762a1bSJed Brown       else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]);
606c4762a1bSJed Brown       if (rdiff > check->tolerance) {
607c4762a1bSJed Brown         tmp        = xa[i];
608c4762a1bSJed Brown         xa[i]      = .5*(xa[i] + xa_last[i]);
609c4762a1bSJed Brown         *changed_x = PETSC_TRUE;
610*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"  Altering entry %D: x=%g, x_last=%g, diff=%g, x_new=%g\n",i,(double)PetscAbsScalar(tmp),(double)PetscAbsScalar(xa_last[i]),(double)rdiff,(double)PetscAbsScalar(xa[i])));
611c4762a1bSJed Brown       }
612c4762a1bSJed Brown     }
613*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,check->last_step,&xa_last));
614*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,x,&xa));
615c4762a1bSJed Brown   }
616*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,check->last_step));
617c4762a1bSJed Brown   PetscFunctionReturn(0);
618c4762a1bSJed Brown }
619c4762a1bSJed Brown 
620c4762a1bSJed Brown /* ------------------------------------------------------------------- */
621c4762a1bSJed Brown /*
622c4762a1bSJed Brown    PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used
623c4762a1bSJed Brown    e.g,
624c4762a1bSJed 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
625c4762a1bSJed Brown    Set by SNESLineSearchSetPostCheck().
626c4762a1bSJed Brown 
627c4762a1bSJed Brown    Input Parameters:
628c4762a1bSJed Brown    linesearch - the LineSearch context
629c4762a1bSJed Brown    xcurrent - current solution
630c4762a1bSJed Brown    y - search direction and length
631c4762a1bSJed Brown    x    - the new candidate iterate
632c4762a1bSJed Brown 
633c4762a1bSJed Brown    Output Parameters:
634c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
635c4762a1bSJed Brown    x    - current iterate (possibly modified)
636c4762a1bSJed Brown 
637c4762a1bSJed Brown  */
638c4762a1bSJed Brown PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
639c4762a1bSJed Brown {
640c4762a1bSJed Brown   SetSubKSPCtx   *check;
641c4762a1bSJed Brown   PetscInt       iter,its,sub_its,maxit;
642c4762a1bSJed Brown   KSP            ksp,sub_ksp,*sub_ksps;
643c4762a1bSJed Brown   PC             pc;
644c4762a1bSJed Brown   PetscReal      ksp_ratio;
645c4762a1bSJed Brown   SNES           snes;
646c4762a1bSJed Brown 
647c4762a1bSJed Brown   PetscFunctionBeginUser;
648*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESLineSearchGetSNES(linesearch, &snes));
649c4762a1bSJed Brown   check   = (SetSubKSPCtx*)ctx;
650*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&iter));
651*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
652*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
653*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCBJacobiGetSubKSP(pc,NULL,NULL,&sub_ksps));
654c4762a1bSJed Brown   sub_ksp = sub_ksps[0];
655*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(ksp,&its));      /* outer KSP iteration number */
656*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(sub_ksp,&sub_its)); /* inner KSP iteration number */
657c4762a1bSJed Brown 
658c4762a1bSJed Brown   if (iter) {
659*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"    ...PostCheck snes iteration %D, ksp_it %D %D, subksp_it %D\n",iter,check->its0,its,sub_its));
660c4762a1bSJed Brown     ksp_ratio = ((PetscReal)(its))/check->its0;
661c4762a1bSJed Brown     maxit     = (PetscInt)(ksp_ratio*sub_its + 0.5);
662c4762a1bSJed Brown     if (maxit < 2) maxit = 2;
663*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit));
664*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"    ...ksp_ratio %g, new maxit %D\n\n",(double)ksp_ratio,maxit));
665c4762a1bSJed Brown   }
666c4762a1bSJed Brown   check->its0 = its; /* save current outer KSP iteration number */
667c4762a1bSJed Brown   PetscFunctionReturn(0);
668c4762a1bSJed Brown }
669c4762a1bSJed Brown 
670c4762a1bSJed Brown /* ------------------------------------------------------------------- */
671c4762a1bSJed Brown /*
672c4762a1bSJed Brown    MatrixFreePreconditioner - This routine demonstrates the use of a
673c4762a1bSJed Brown    user-provided preconditioner.  This code implements just the null
674c4762a1bSJed Brown    preconditioner, which of course is not recommended for general use.
675c4762a1bSJed Brown 
676c4762a1bSJed Brown    Input Parameters:
677c4762a1bSJed Brown +  pc - preconditioner
678c4762a1bSJed Brown -  x - input vector
679c4762a1bSJed Brown 
680c4762a1bSJed Brown    Output Parameter:
681c4762a1bSJed Brown .  y - preconditioned vector
682c4762a1bSJed Brown */
683c4762a1bSJed Brown PetscErrorCode MatrixFreePreconditioner(PC pc,Vec x,Vec y)
684c4762a1bSJed Brown {
685*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,y));
686c4762a1bSJed Brown   return 0;
687c4762a1bSJed Brown }
688c4762a1bSJed Brown 
689c4762a1bSJed Brown /*TEST
690c4762a1bSJed Brown 
691c4762a1bSJed Brown    test:
692c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
693c4762a1bSJed Brown 
694c4762a1bSJed Brown    test:
695c4762a1bSJed Brown       suffix: 2
696c4762a1bSJed Brown       nsize: 3
697c4762a1bSJed Brown       args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
698c4762a1bSJed Brown 
699c4762a1bSJed Brown    test:
700c4762a1bSJed Brown       suffix: 3
701c4762a1bSJed Brown       nsize: 2
702c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
703c4762a1bSJed Brown 
704c4762a1bSJed Brown    test:
705c4762a1bSJed Brown       suffix: 4
706c4762a1bSJed Brown       args: -nox -pre_check_iterates -post_check_iterates
707c4762a1bSJed Brown 
708c4762a1bSJed Brown    test:
709c4762a1bSJed Brown       suffix: 5
710c4762a1bSJed Brown       requires: double !complex !single
711c4762a1bSJed Brown       nsize: 2
712c4762a1bSJed Brown       args: -nox -snes_test_jacobian  -snes_test_jacobian_view
713c4762a1bSJed Brown 
714c4762a1bSJed Brown    test:
715c4762a1bSJed Brown       suffix: 6
716c4762a1bSJed Brown       requires: double !complex !single
717c4762a1bSJed Brown       nsize: 4
718c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1
719c4762a1bSJed Brown 
720c4762a1bSJed Brown    test:
721c4762a1bSJed Brown       suffix: 7
722c4762a1bSJed Brown       requires: double !complex !single
723c4762a1bSJed Brown       nsize: 4
724c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1
725c4762a1bSJed Brown 
726c4762a1bSJed Brown    test:
727c4762a1bSJed Brown       suffix: 8
728c4762a1bSJed Brown       requires: double !complex !single
729c4762a1bSJed Brown       nsize: 4
730c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1
731c4762a1bSJed Brown 
732c4762a1bSJed Brown    test:
733c4762a1bSJed Brown       suffix: 9
734c4762a1bSJed Brown       requires: double !complex !single
735c4762a1bSJed Brown       nsize: 4
736c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1
737c4762a1bSJed Brown 
738c4762a1bSJed Brown    test:
739c4762a1bSJed Brown       suffix: 10
740c4762a1bSJed Brown       requires: double !complex !single
741c4762a1bSJed Brown       nsize: 4
742c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1
743c4762a1bSJed Brown 
744c4762a1bSJed Brown    test:
745c4762a1bSJed Brown       suffix: 11
746c4762a1bSJed Brown       requires: double !complex !single
747c4762a1bSJed Brown       nsize: 4
74857715debSLisandro 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
749c4762a1bSJed Brown 
750c0558f20SBarry Smith    test:
751c0558f20SBarry Smith       suffix: 12
752c0558f20SBarry Smith       args: -view_initial
753c0558f20SBarry Smith       filter: grep -v "type:"
754c0558f20SBarry Smith 
75541ba4c6cSHeeho Park    test:
75641ba4c6cSHeeho Park       suffix: 13
75741ba4c6cSHeeho Park       requires: double !complex !single
75841ba4c6cSHeeho Park       nsize: 4
75941ba4c6cSHeeho Park       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1
76041ba4c6cSHeeho Park 
761c4762a1bSJed Brown TEST*/
762