xref: /petsc/src/snes/tutorials/ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscInt       its,N = 5,i,maxit,maxf,xs,xm;
93c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
94c0558f20SBarry Smith   PetscBool      flg,viewinitial = PETSC_FALSE;
95c4762a1bSJed Brown 
96*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
975f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank));
985f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
100c4762a1bSJed Brown   ctx.h = 1.0/(N-1);
101c4762a1bSJed Brown   ctx.sjerr = PETSC_FALSE;
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_jacobian_domain_error",&ctx.sjerr,NULL));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown      Create nonlinear solver context
107c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108c4762a1bSJed Brown 
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
113c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*
116c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
117c4762a1bSJed Brown   */
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(ctx.da));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(ctx.da));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /*
123c4762a1bSJed Brown      Extract global and local vectors from DMDA; then duplicate for remaining
124c4762a1bSJed Brown      vectors that are the same types
125c4762a1bSJed Brown   */
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(ctx.da,&x));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F)); ctx.F = F;
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /*
132c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
133c4762a1bSJed Brown      solver needs to compute the nonlinear function, it will call this
134c4762a1bSJed Brown      routine.
135c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
136c4762a1bSJed Brown         context that provides application-specific data for the
137c4762a1bSJed Brown         function evaluation routine.
138c4762a1bSJed Brown   */
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormFunction,&ctx));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
143c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144c4762a1bSJed Brown 
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(J,3,NULL,3,NULL));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /*
152c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
153c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
154c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
155c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
156c4762a1bSJed Brown         context that provides application-specific data for the
157c4762a1bSJed Brown         Jacobian evaluation routine.
158c4762a1bSJed Brown   */
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,&ctx));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Optionally allow user-provided preconditioner
163c4762a1bSJed Brown    */
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-user_precond",&flg));
165c4762a1bSJed Brown   if (flg) {
166c4762a1bSJed Brown     KSP ksp;
167c4762a1bSJed Brown     PC  pc;
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetKSP(snes,&ksp));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCSHELL));
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellSetApply(pc,MatrixFreePreconditioner));
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
176c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /*
179c4762a1bSJed Brown      Set an optional user-defined monitoring routine
180c4762a1bSJed Brown   */
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorSet(snes,Monitor,&monP,0));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /*
185c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
186c4762a1bSJed Brown   */
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution"));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /*
191c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
192c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
193c4762a1bSJed Brown   */
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /*
197c4762a1bSJed Brown      Set an optional user-defined routine to check the validity of candidate
198c4762a1bSJed Brown      iterates that are determined by line search methods
199c4762a1bSJed Brown   */
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLineSearch(snes, &linesearch));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-post_check_iterates",&post_check));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   if (post_check) {
2045f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating post step checking routine\n"));
2055f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESLineSearchSetPostCheck(linesearch,PostCheck,&checkP));
2065f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&(checkP.last_step)));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown     checkP.tolerance = 1.0;
209c4762a1bSJed Brown     checkP.user      = &ctx;
210c4762a1bSJed Brown 
2115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-check_tol",&checkP.tolerance,NULL));
212c4762a1bSJed Brown   }
213c4762a1bSJed Brown 
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-post_setsubksp",&post_setsubksp));
215c4762a1bSJed Brown   if (post_setsubksp) {
2165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating post setsubksp\n"));
2175f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESLineSearchSetPostCheck(linesearch,PostSetSubKSP,&checkP1));
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown 
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-pre_check_iterates",&pre_check));
221c4762a1bSJed Brown   if (pre_check) {
2225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Activating pre step checking routine\n"));
2235f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESLineSearchSetPreCheck(linesearch,PreCheck,&checkP));
224c4762a1bSJed Brown   }
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /*
227c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
228c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
229c4762a1bSJed Brown      the option -snes_view
230c4762a1bSJed Brown   */
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
2325f80ce2aSJacob 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));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235c4762a1bSJed Brown      Initialize application:
236c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
237c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /*
240c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
241c4762a1bSJed Brown        xs, xm - starting grid index, width of local grid (no ghost points)
242c4762a1bSJed Brown   */
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(ctx.da,&xs,NULL,NULL,&xm,NULL,NULL));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /*
246c4762a1bSJed Brown      Get pointers to vector data
247c4762a1bSJed Brown   */
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(ctx.da,F,&FF));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(ctx.da,U,&UU));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /*
252c4762a1bSJed Brown      Compute local vector entries
253c4762a1bSJed Brown   */
254c4762a1bSJed Brown   xp = ctx.h*xs;
255c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
256c4762a1bSJed Brown     FF[i] = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
257c4762a1bSJed Brown     UU[i] = xp*xp*xp;
258c4762a1bSJed Brown     xp   += ctx.h;
259c4762a1bSJed Brown   }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /*
262c4762a1bSJed Brown      Restore vectors
263c4762a1bSJed Brown   */
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(ctx.da,F,&FF));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(ctx.da,U,&UU));
266c0558f20SBarry Smith   if (viewinitial) {
2675f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(U,0));
2685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(F,0));
269c0558f20SBarry Smith   }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
273c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /*
276c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
277c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
278c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
279c4762a1bSJed Brown      this vector to zero by calling VecSet().
280c4762a1bSJed Brown   */
2815f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287c4762a1bSJed Brown      Check solution and clean up
288c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   /*
291c4762a1bSJed Brown      Check the error
292c4762a1bSJed Brown   */
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,none,U));
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its));
296c4762a1bSJed Brown   if (ctx.sjerr) {
297c4762a1bSJed Brown     SNESType       snestype;
2985f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetType(snes,&snestype));
2995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"SNES Type %s\n",snestype));
300c4762a1bSJed Brown   }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /*
303c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
304c4762a1bSJed Brown      are no longer needed.
305c4762a1bSJed Brown   */
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&monP.viewer));
3075f80ce2aSJacob Faibussowitsch   if (post_check) CHKERRQ(VecDestroy(&checkP.last_step));
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&F));
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&ctx.da));
315*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
316*b122ec5aSJacob Faibussowitsch   return 0;
317c4762a1bSJed Brown }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown /* ------------------------------------------------------------------- */
320c4762a1bSJed Brown /*
321c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
322c4762a1bSJed Brown 
323c4762a1bSJed Brown    Input/Output Parameter:
324c4762a1bSJed Brown .  x - the solution vector
325c4762a1bSJed Brown */
326c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
327c4762a1bSJed Brown {
328c4762a1bSJed Brown   PetscScalar    pfive = .50;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBeginUser;
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
332c4762a1bSJed Brown   PetscFunctionReturn(0);
333c4762a1bSJed Brown }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown /* ------------------------------------------------------------------- */
336c4762a1bSJed Brown /*
337c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
338c4762a1bSJed Brown 
339c4762a1bSJed Brown    Input Parameters:
340c4762a1bSJed Brown .  snes - the SNES context
341c4762a1bSJed Brown .  x - input vector
342c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
343c4762a1bSJed Brown 
344c4762a1bSJed Brown    Output Parameter:
345c4762a1bSJed Brown .  f - function vector
346c4762a1bSJed Brown 
347c4762a1bSJed Brown    Note:
348c4762a1bSJed Brown    The user-defined context can contain any application-specific
349c4762a1bSJed Brown    data needed for the function evaluation.
350c4762a1bSJed Brown */
351c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
352c4762a1bSJed Brown {
353c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx*) ctx;
354c4762a1bSJed Brown   DM             da    = user->da;
355c4762a1bSJed Brown   PetscScalar    *xx,*ff,*FF,d;
356c4762a1bSJed Brown   PetscInt       i,M,xs,xm;
357c4762a1bSJed Brown   Vec            xlocal;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   PetscFunctionBeginUser;
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&xlocal));
361c4762a1bSJed Brown   /*
362c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
363c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
364c4762a1bSJed Brown      By placing code between these two statements, computations can
365c4762a1bSJed Brown      be done while messages are in transition.
366c4762a1bSJed Brown   */
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   /*
371c4762a1bSJed Brown      Get pointers to vector data.
372c4762a1bSJed Brown        - The vector xlocal includes ghost point; the vectors x and f do
373c4762a1bSJed Brown          NOT include ghost points.
374c4762a1bSJed Brown        - Using DMDAVecGetArray() allows accessing the values using global ordering
375c4762a1bSJed Brown   */
3765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,xlocal,&xx));
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,f,&ff));
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,user->F,&FF));
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   /*
381c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
382c4762a1bSJed Brown        xs, xm  - starting grid index, width of local grid (no ghost points)
383c4762a1bSJed Brown   */
3845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
3855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /*
388c4762a1bSJed Brown      Set function values for boundary points; define local interior grid point range:
389c4762a1bSJed Brown         xsi - starting interior grid index
390c4762a1bSJed Brown         xei - ending interior grid index
391c4762a1bSJed Brown   */
392c4762a1bSJed Brown   if (xs == 0) { /* left boundary */
393c4762a1bSJed Brown     ff[0] = xx[0];
394c4762a1bSJed Brown     xs++;xm--;
395c4762a1bSJed Brown   }
396c4762a1bSJed Brown   if (xs+xm == M) {  /* right boundary */
397c4762a1bSJed Brown     ff[xs+xm-1] = xx[xs+xm-1] - 1.0;
398c4762a1bSJed Brown     xm--;
399c4762a1bSJed Brown   }
400c4762a1bSJed Brown 
401c4762a1bSJed Brown   /*
402c4762a1bSJed Brown      Compute function over locally owned part of the grid (interior points only)
403c4762a1bSJed Brown   */
404c4762a1bSJed Brown   d = 1.0/(user->h*user->h);
405c4762a1bSJed 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];
406c4762a1bSJed Brown 
407c4762a1bSJed Brown   /*
408c4762a1bSJed Brown      Restore vectors
409c4762a1bSJed Brown   */
4105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,xlocal,&xx));
4115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,f,&ff));
4125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,user->F,&FF));
4135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&xlocal));
414c4762a1bSJed Brown   PetscFunctionReturn(0);
415c4762a1bSJed Brown }
416c4762a1bSJed Brown 
417c4762a1bSJed Brown /* ------------------------------------------------------------------- */
418c4762a1bSJed Brown /*
419c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
420c4762a1bSJed Brown 
421c4762a1bSJed Brown    Input Parameters:
422c4762a1bSJed Brown .  snes - the SNES context
423c4762a1bSJed Brown .  x - input vector
424c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
425c4762a1bSJed Brown 
426c4762a1bSJed Brown    Output Parameters:
427c4762a1bSJed Brown .  jac - Jacobian matrix
428c4762a1bSJed Brown .  B - optionally different preconditioning matrix
429c4762a1bSJed Brown .  flag - flag indicating matrix structure
430c4762a1bSJed Brown */
431c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
432c4762a1bSJed Brown {
433c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx*) ctx;
434c4762a1bSJed Brown   PetscScalar    *xx,d,A[3];
435c4762a1bSJed Brown   PetscInt       i,j[3],M,xs,xm;
436c4762a1bSJed Brown   DM             da = user->da;
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   PetscFunctionBeginUser;
439c4762a1bSJed Brown   if (user->sjerr) {
4405f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetJacobianDomainError(snes));
441c4762a1bSJed Brown     PetscFunctionReturn(0);
442c4762a1bSJed Brown   }
443c4762a1bSJed Brown   /*
444c4762a1bSJed Brown      Get pointer to vector data
445c4762a1bSJed Brown   */
4465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,x,&xx));
4475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /*
450c4762a1bSJed Brown     Get range of locally owned matrix
451c4762a1bSJed Brown   */
4525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   /*
455c4762a1bSJed Brown      Determine starting and ending local indices for interior grid points.
456c4762a1bSJed Brown      Set Jacobian entries for boundary points.
457c4762a1bSJed Brown   */
458c4762a1bSJed Brown 
459c4762a1bSJed Brown   if (xs == 0) {  /* left boundary */
460c4762a1bSJed Brown     i = 0; A[0] = 1.0;
461c4762a1bSJed Brown 
4625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
463c4762a1bSJed Brown     xs++;xm--;
464c4762a1bSJed Brown   }
465c4762a1bSJed Brown   if (xs+xm == M) { /* right boundary */
466c4762a1bSJed Brown     i    = M-1;
467c4762a1bSJed Brown     A[0] = 1.0;
4685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
469c4762a1bSJed Brown     xm--;
470c4762a1bSJed Brown   }
471c4762a1bSJed Brown 
472c4762a1bSJed Brown   /*
473c4762a1bSJed Brown      Interior grid points
474c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
475c4762a1bSJed Brown         row at once.
476c4762a1bSJed Brown   */
477c4762a1bSJed Brown   d = 1.0/(user->h*user->h);
478c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
479c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
480c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
4815f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
482c4762a1bSJed Brown   }
483c4762a1bSJed Brown 
484c4762a1bSJed Brown   /*
485c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
486c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
487c4762a1bSJed Brown      By placing code between these two statements, computations can be
488c4762a1bSJed Brown      done while messages are in transition.
489c4762a1bSJed Brown 
490c4762a1bSJed Brown      Also, restore vector.
491c4762a1bSJed Brown   */
492c4762a1bSJed Brown 
4935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
4945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,x,&xx));
4955f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   PetscFunctionReturn(0);
498c4762a1bSJed Brown }
499c4762a1bSJed Brown 
500c4762a1bSJed Brown /* ------------------------------------------------------------------- */
501c4762a1bSJed Brown /*
502c4762a1bSJed Brown    Monitor - Optional user-defined monitoring routine that views the
503c4762a1bSJed Brown    current iterate with an x-window plot. Set by SNESMonitorSet().
504c4762a1bSJed Brown 
505c4762a1bSJed Brown    Input Parameters:
506c4762a1bSJed Brown    snes - the SNES context
507c4762a1bSJed Brown    its - iteration number
508c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
509c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
510c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
511c4762a1bSJed Brown 
512c4762a1bSJed Brown    Note:
513c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
514c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
515c4762a1bSJed Brown  */
516c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
517c4762a1bSJed Brown {
518c4762a1bSJed Brown   MonitorCtx     *monP = (MonitorCtx*) ctx;
519c4762a1bSJed Brown   Vec            x;
520c4762a1bSJed Brown 
521c4762a1bSJed Brown   PetscFunctionBeginUser;
5225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,(double)fnorm));
5235f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes,&x));
5245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,monP->viewer));
525c4762a1bSJed Brown   PetscFunctionReturn(0);
526c4762a1bSJed Brown }
527c4762a1bSJed Brown 
528c4762a1bSJed Brown /* ------------------------------------------------------------------- */
529c4762a1bSJed Brown /*
530c4762a1bSJed Brown    PreCheck - Optional user-defined routine that checks the validity of
531c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPreCheck().
532c4762a1bSJed Brown 
533c4762a1bSJed Brown    Input Parameters:
534c4762a1bSJed Brown    snes - the SNES context
535c4762a1bSJed Brown    xcurrent - current solution
536c4762a1bSJed Brown    y - search direction and length
537c4762a1bSJed Brown 
538c4762a1bSJed Brown    Output Parameters:
539c4762a1bSJed Brown    y         - proposed step (search direction and length) (possibly changed)
540c4762a1bSJed Brown    changed_y - tells if the step has changed or not
541c4762a1bSJed Brown  */
542c4762a1bSJed Brown PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx)
543c4762a1bSJed Brown {
544c4762a1bSJed Brown   PetscFunctionBeginUser;
545c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
546c4762a1bSJed Brown   PetscFunctionReturn(0);
547c4762a1bSJed Brown }
548c4762a1bSJed Brown 
549c4762a1bSJed Brown /* ------------------------------------------------------------------- */
550c4762a1bSJed Brown /*
551c4762a1bSJed Brown    PostCheck - Optional user-defined routine that checks the validity of
552c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPostCheck().
553c4762a1bSJed Brown 
554c4762a1bSJed Brown    Input Parameters:
555c4762a1bSJed Brown    snes - the SNES context
556c4762a1bSJed Brown    ctx  - optional user-defined context for private data for the
557c4762a1bSJed Brown           monitor routine, as set by SNESLineSearchSetPostCheck()
558c4762a1bSJed Brown    xcurrent - current solution
559c4762a1bSJed Brown    y - search direction and length
560c4762a1bSJed Brown    x    - the new candidate iterate
561c4762a1bSJed Brown 
562c4762a1bSJed Brown    Output Parameters:
563c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
564c4762a1bSJed Brown    x    - current iterate (possibly modified)
565c4762a1bSJed Brown 
566c4762a1bSJed Brown  */
567c4762a1bSJed Brown PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
568c4762a1bSJed Brown {
569c4762a1bSJed Brown   PetscInt       i,iter,xs,xm;
570c4762a1bSJed Brown   StepCheckCtx   *check;
571c4762a1bSJed Brown   ApplicationCtx *user;
572c4762a1bSJed Brown   PetscScalar    *xa,*xa_last,tmp;
573c4762a1bSJed Brown   PetscReal      rdiff;
574c4762a1bSJed Brown   DM             da;
575c4762a1bSJed Brown   SNES           snes;
576c4762a1bSJed Brown 
577c4762a1bSJed Brown   PetscFunctionBeginUser;
578c4762a1bSJed Brown   *changed_x = PETSC_FALSE;
579c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
580c4762a1bSJed Brown 
5815f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESLineSearchGetSNES(linesearch, &snes));
582c4762a1bSJed Brown   check = (StepCheckCtx*)ctx;
583c4762a1bSJed Brown   user  = check->user;
5845f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&iter));
585c4762a1bSJed Brown 
586c4762a1bSJed Brown   /* iteration 1 indicates we are working on the second iteration */
587c4762a1bSJed Brown   if (iter > 0) {
588c4762a1bSJed Brown     da   = user->da;
5895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",iter,(double)check->tolerance));
590c4762a1bSJed Brown 
591c4762a1bSJed Brown     /* Access local array data */
5925f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,check->last_step,&xa_last));
5935f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,x,&xa));
5945f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
595c4762a1bSJed Brown 
596c4762a1bSJed Brown     /*
597c4762a1bSJed Brown        If we fail the user-defined check for validity of the candidate iterate,
598c4762a1bSJed Brown        then modify the iterate as we like.  (Note that the particular modification
599c4762a1bSJed Brown        below is intended simply to demonstrate how to manipulate this data, not
600c4762a1bSJed Brown        as a meaningful or appropriate choice.)
601c4762a1bSJed Brown     */
602c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
603c4762a1bSJed Brown       if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance;
604c4762a1bSJed Brown       else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]);
605c4762a1bSJed Brown       if (rdiff > check->tolerance) {
606c4762a1bSJed Brown         tmp        = xa[i];
607c4762a1bSJed Brown         xa[i]      = .5*(xa[i] + xa_last[i]);
608c4762a1bSJed Brown         *changed_x = PETSC_TRUE;
6095f80ce2aSJacob 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])));
610c4762a1bSJed Brown       }
611c4762a1bSJed Brown     }
6125f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,check->last_step,&xa_last));
6135f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,x,&xa));
614c4762a1bSJed Brown   }
6155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,check->last_step));
616c4762a1bSJed Brown   PetscFunctionReturn(0);
617c4762a1bSJed Brown }
618c4762a1bSJed Brown 
619c4762a1bSJed Brown /* ------------------------------------------------------------------- */
620c4762a1bSJed Brown /*
621c4762a1bSJed Brown    PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used
622c4762a1bSJed Brown    e.g,
623c4762a1bSJed 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
624c4762a1bSJed Brown    Set by SNESLineSearchSetPostCheck().
625c4762a1bSJed Brown 
626c4762a1bSJed Brown    Input Parameters:
627c4762a1bSJed Brown    linesearch - the LineSearch context
628c4762a1bSJed Brown    xcurrent - current solution
629c4762a1bSJed Brown    y - search direction and length
630c4762a1bSJed Brown    x    - the new candidate iterate
631c4762a1bSJed Brown 
632c4762a1bSJed Brown    Output Parameters:
633c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
634c4762a1bSJed Brown    x    - current iterate (possibly modified)
635c4762a1bSJed Brown 
636c4762a1bSJed Brown  */
637c4762a1bSJed Brown PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
638c4762a1bSJed Brown {
639c4762a1bSJed Brown   SetSubKSPCtx   *check;
640c4762a1bSJed Brown   PetscInt       iter,its,sub_its,maxit;
641c4762a1bSJed Brown   KSP            ksp,sub_ksp,*sub_ksps;
642c4762a1bSJed Brown   PC             pc;
643c4762a1bSJed Brown   PetscReal      ksp_ratio;
644c4762a1bSJed Brown   SNES           snes;
645c4762a1bSJed Brown 
646c4762a1bSJed Brown   PetscFunctionBeginUser;
6475f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESLineSearchGetSNES(linesearch, &snes));
648c4762a1bSJed Brown   check   = (SetSubKSPCtx*)ctx;
6495f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&iter));
6505f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
6515f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
6525f80ce2aSJacob Faibussowitsch   CHKERRQ(PCBJacobiGetSubKSP(pc,NULL,NULL,&sub_ksps));
653c4762a1bSJed Brown   sub_ksp = sub_ksps[0];
6545f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(ksp,&its));      /* outer KSP iteration number */
6555f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(sub_ksp,&sub_its)); /* inner KSP iteration number */
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   if (iter) {
6585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"    ...PostCheck snes iteration %D, ksp_it %D %D, subksp_it %D\n",iter,check->its0,its,sub_its));
659c4762a1bSJed Brown     ksp_ratio = ((PetscReal)(its))/check->its0;
660c4762a1bSJed Brown     maxit     = (PetscInt)(ksp_ratio*sub_its + 0.5);
661c4762a1bSJed Brown     if (maxit < 2) maxit = 2;
6625f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit));
6635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"    ...ksp_ratio %g, new maxit %D\n\n",(double)ksp_ratio,maxit));
664c4762a1bSJed Brown   }
665c4762a1bSJed Brown   check->its0 = its; /* save current outer KSP iteration number */
666c4762a1bSJed Brown   PetscFunctionReturn(0);
667c4762a1bSJed Brown }
668c4762a1bSJed Brown 
669c4762a1bSJed Brown /* ------------------------------------------------------------------- */
670c4762a1bSJed Brown /*
671c4762a1bSJed Brown    MatrixFreePreconditioner - This routine demonstrates the use of a
672c4762a1bSJed Brown    user-provided preconditioner.  This code implements just the null
673c4762a1bSJed Brown    preconditioner, which of course is not recommended for general use.
674c4762a1bSJed Brown 
675c4762a1bSJed Brown    Input Parameters:
676c4762a1bSJed Brown +  pc - preconditioner
677c4762a1bSJed Brown -  x - input vector
678c4762a1bSJed Brown 
679c4762a1bSJed Brown    Output Parameter:
680c4762a1bSJed Brown .  y - preconditioned vector
681c4762a1bSJed Brown */
682c4762a1bSJed Brown PetscErrorCode MatrixFreePreconditioner(PC pc,Vec x,Vec y)
683c4762a1bSJed Brown {
6845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,y));
685c4762a1bSJed Brown   return 0;
686c4762a1bSJed Brown }
687c4762a1bSJed Brown 
688c4762a1bSJed Brown /*TEST
689c4762a1bSJed Brown 
690c4762a1bSJed Brown    test:
691c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
692c4762a1bSJed Brown 
693c4762a1bSJed Brown    test:
694c4762a1bSJed Brown       suffix: 2
695c4762a1bSJed Brown       nsize: 3
696c4762a1bSJed Brown       args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
697c4762a1bSJed Brown 
698c4762a1bSJed Brown    test:
699c4762a1bSJed Brown       suffix: 3
700c4762a1bSJed Brown       nsize: 2
701c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
702c4762a1bSJed Brown 
703c4762a1bSJed Brown    test:
704c4762a1bSJed Brown       suffix: 4
705c4762a1bSJed Brown       args: -nox -pre_check_iterates -post_check_iterates
706c4762a1bSJed Brown 
707c4762a1bSJed Brown    test:
708c4762a1bSJed Brown       suffix: 5
709c4762a1bSJed Brown       requires: double !complex !single
710c4762a1bSJed Brown       nsize: 2
711c4762a1bSJed Brown       args: -nox -snes_test_jacobian  -snes_test_jacobian_view
712c4762a1bSJed Brown 
713c4762a1bSJed Brown    test:
714c4762a1bSJed Brown       suffix: 6
715c4762a1bSJed Brown       requires: double !complex !single
716c4762a1bSJed Brown       nsize: 4
717c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1
718c4762a1bSJed Brown 
719c4762a1bSJed Brown    test:
720c4762a1bSJed Brown       suffix: 7
721c4762a1bSJed Brown       requires: double !complex !single
722c4762a1bSJed Brown       nsize: 4
723c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1
724c4762a1bSJed Brown 
725c4762a1bSJed Brown    test:
726c4762a1bSJed Brown       suffix: 8
727c4762a1bSJed Brown       requires: double !complex !single
728c4762a1bSJed Brown       nsize: 4
729c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1
730c4762a1bSJed Brown 
731c4762a1bSJed Brown    test:
732c4762a1bSJed Brown       suffix: 9
733c4762a1bSJed Brown       requires: double !complex !single
734c4762a1bSJed Brown       nsize: 4
735c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1
736c4762a1bSJed Brown 
737c4762a1bSJed Brown    test:
738c4762a1bSJed Brown       suffix: 10
739c4762a1bSJed Brown       requires: double !complex !single
740c4762a1bSJed Brown       nsize: 4
741c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1
742c4762a1bSJed Brown 
743c4762a1bSJed Brown    test:
744c4762a1bSJed Brown       suffix: 11
745c4762a1bSJed Brown       requires: double !complex !single
746c4762a1bSJed Brown       nsize: 4
74757715debSLisandro 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
748c4762a1bSJed Brown 
749c0558f20SBarry Smith    test:
750c0558f20SBarry Smith       suffix: 12
751c0558f20SBarry Smith       args: -view_initial
752c0558f20SBarry Smith       filter: grep -v "type:"
753c0558f20SBarry Smith 
75441ba4c6cSHeeho Park    test:
75541ba4c6cSHeeho Park       suffix: 13
75641ba4c6cSHeeho Park       requires: double !complex !single
75741ba4c6cSHeeho Park       nsize: 4
75841ba4c6cSHeeho Park       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1
75941ba4c6cSHeeho Park 
760c4762a1bSJed Brown TEST*/
761