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