xref: /petsc/src/snes/tutorials/ex3.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 */
449371c9d4SSatish Balay typedef struct {
45c4762a1bSJed Brown   DM          da;    /* distributed array */
46dd8e379bSPierre Jolivet   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 */
569371c9d4SSatish Balay 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 */
649371c9d4SSatish Balay 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 
709371c9d4SSatish Balay typedef struct {
71da81f932SPierre Jolivet   PetscInt its0; /* num of previous outer KSP iterations */
72c4762a1bSJed Brown } SetSubKSPCtx;
73c4762a1bSJed Brown 
main(int argc,char ** argv)74d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
75d71ae5a4SJacob Faibussowitsch {
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 
90327415f7SBarry Smith   PetscFunctionBeginUser;
91c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &ctx.rank));
939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &ctx.size));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL));
95c4762a1bSJed Brown   ctx.h     = 1.0 / (N - 1);
96c4762a1bSJed Brown   ctx.sjerr = PETSC_FALSE;
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_jacobian_domain_error", &ctx.sjerr, NULL));
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown      Create nonlinear solver context
102c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
108c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /*
111c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
112c4762a1bSJed Brown   */
1139566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, N, 1, 1, NULL, &ctx.da));
1149566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(ctx.da));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetUp(ctx.da));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /*
118c4762a1bSJed Brown      Extract global and local vectors from DMDA; then duplicate for remaining
119c4762a1bSJed Brown      vectors that are the same types
120c4762a1bSJed Brown   */
1219566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(ctx.da, &x));
1229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
12311486bccSBarry Smith   PetscCall(VecDuplicate(x, &F));
12411486bccSBarry Smith   ctx.F = F;
1259566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /*
128c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
129c4762a1bSJed Brown      solver needs to compute the nonlinear function, it will call this
130c4762a1bSJed Brown      routine.
131c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
132c4762a1bSJed Brown         context that provides application-specific data for the
133c4762a1bSJed Brown         function evaluation routine.
134c4762a1bSJed Brown   */
1359566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, &ctx));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
139c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
1429566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, N, N));
1439566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
1449566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
1459566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(J, 3, NULL, 3, NULL));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
149c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
150c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
151c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
152c4762a1bSJed Brown         context that provides application-specific data for the
153c4762a1bSJed Brown         Jacobian evaluation routine.
154c4762a1bSJed Brown   */
1559566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &ctx));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /*
158c4762a1bSJed Brown      Optionally allow user-provided preconditioner
159c4762a1bSJed Brown    */
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-user_precond", &flg));
161c4762a1bSJed Brown   if (flg) {
162c4762a1bSJed Brown     KSP ksp;
163c4762a1bSJed Brown     PC  pc;
1649566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
1659566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
1669566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCSHELL));
1679566063dSJacob Faibussowitsch     PetscCall(PCShellSetApply(pc, MatrixFreePreconditioner));
168c4762a1bSJed Brown   }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
172c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /*
175c4762a1bSJed Brown      Set an optional user-defined monitoring routine
176c4762a1bSJed Brown   */
1779566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, 0, 0, 0, 400, 400, &monP.viewer));
1789566063dSJacob Faibussowitsch   PetscCall(SNESMonitorSet(snes, Monitor, &monP, 0));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /*
181c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
182c4762a1bSJed Brown   */
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution"));
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution"));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /*
187c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
188c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
189c4762a1bSJed Brown   */
1909566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /*
193c4762a1bSJed Brown      Set an optional user-defined routine to check the validity of candidate
194c4762a1bSJed Brown      iterates that are determined by line search methods
195c4762a1bSJed Brown   */
1969566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
1979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-post_check_iterates", &post_check));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   if (post_check) {
2009566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating post step checking routine\n"));
2019566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPostCheck(linesearch, PostCheck, &checkP));
202f4f49eeaSPierre Jolivet     PetscCall(VecDuplicate(x, &checkP.last_step));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown     checkP.tolerance = 1.0;
205c4762a1bSJed Brown     checkP.user      = &ctx;
206c4762a1bSJed Brown 
2079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-check_tol", &checkP.tolerance, NULL));
208c4762a1bSJed Brown   }
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-post_setsubksp", &post_setsubksp));
211c4762a1bSJed Brown   if (post_setsubksp) {
2129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating post setsubksp\n"));
2139566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPostCheck(linesearch, PostSetSubKSP, &checkP1));
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown 
2169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-pre_check_iterates", &pre_check));
217c4762a1bSJed Brown   if (pre_check) {
2189566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating pre step checking routine\n"));
2199566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetPreCheck(linesearch, PreCheck, &checkP));
220c4762a1bSJed Brown   }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /*
223c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
224c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
225c4762a1bSJed Brown      the option -snes_view
226c4762a1bSJed Brown   */
2279566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf));
22863a3b9bcSJacob 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));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231c4762a1bSJed Brown      Initialize application:
232dd8e379bSPierre Jolivet      Store right-hand side of PDE and exact solution
233c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /*
236c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
237c4762a1bSJed Brown        xs, xm - starting grid index, width of local grid (no ghost points)
238c4762a1bSJed Brown   */
2399566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(ctx.da, &xs, NULL, NULL, &xm, NULL, NULL));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /*
242c4762a1bSJed Brown      Get pointers to vector data
243c4762a1bSJed Brown   */
2449566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(ctx.da, F, &FF));
2459566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(ctx.da, U, &UU));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /*
248c4762a1bSJed Brown      Compute local vector entries
249c4762a1bSJed Brown   */
250c4762a1bSJed Brown   xp = ctx.h * xs;
251c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
252c4762a1bSJed Brown     FF[i] = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
253c4762a1bSJed Brown     UU[i] = xp * xp * xp;
254c4762a1bSJed Brown     xp += ctx.h;
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /*
258c4762a1bSJed Brown      Restore vectors
259c4762a1bSJed Brown   */
2609566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(ctx.da, F, &FF));
2619566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(ctx.da, U, &UU));
262c0558f20SBarry Smith   if (viewinitial) {
2639566063dSJacob Faibussowitsch     PetscCall(VecView(U, 0));
2649566063dSJacob Faibussowitsch     PetscCall(VecView(F, 0));
265c0558f20SBarry Smith   }
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
269c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /*
272c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
273c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
274c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
275c4762a1bSJed Brown      this vector to zero by calling VecSet().
276c4762a1bSJed Brown   */
2779566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
2789566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
2799566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
28063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283c4762a1bSJed Brown      Check solution and clean up
284c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /*
287c4762a1bSJed Brown      Check the error
288c4762a1bSJed Brown   */
2899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, none, U));
2909566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
29163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Iterations %" PetscInt_FMT "\n", (double)norm, its));
292c4762a1bSJed Brown   if (ctx.sjerr) {
293c4762a1bSJed Brown     SNESType snestype;
2949566063dSJacob Faibussowitsch     PetscCall(SNESGetType(snes, &snestype));
2959566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES Type %s\n", snestype));
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   /*
299c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
300c4762a1bSJed Brown      are no longer needed.
301c4762a1bSJed Brown   */
3029566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&monP.viewer));
3039566063dSJacob Faibussowitsch   if (post_check) PetscCall(VecDestroy(&checkP.last_step));
3049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F));
3089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
3099566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
3109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx.da));
3119566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
312b122ec5aSJacob Faibussowitsch   return 0;
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown /* ------------------------------------------------------------------- */
316c4762a1bSJed Brown /*
317c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
318c4762a1bSJed Brown 
319c4762a1bSJed Brown    Input/Output Parameter:
320c4762a1bSJed Brown .  x - the solution vector
321c4762a1bSJed Brown */
FormInitialGuess(Vec x)322d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
323d71ae5a4SJacob Faibussowitsch {
324c4762a1bSJed Brown   PetscScalar pfive = .50;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   PetscFunctionBeginUser;
3279566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
3283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
329c4762a1bSJed Brown }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown /* ------------------------------------------------------------------- */
332c4762a1bSJed Brown /*
333c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
334c4762a1bSJed Brown 
335c4762a1bSJed Brown    Input Parameters:
336c4762a1bSJed Brown .  snes - the SNES context
337c4762a1bSJed Brown .  x - input vector
338c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    Output Parameter:
341c4762a1bSJed Brown .  f - function vector
342c4762a1bSJed Brown 
343c4762a1bSJed Brown    Note:
344c4762a1bSJed Brown    The user-defined context can contain any application-specific
345c4762a1bSJed Brown    data needed for the function evaluation.
346c4762a1bSJed Brown */
FormFunction(SNES snes,Vec x,Vec f,PetscCtx ctx)347*2a8381b2SBarry Smith PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, PetscCtx ctx)
348d71ae5a4SJacob Faibussowitsch {
349c4762a1bSJed Brown   ApplicationCtx    *user = (ApplicationCtx *)ctx;
350c4762a1bSJed Brown   DM                 da   = user->da;
35106cf5b03SBarry Smith   PetscScalar       *ff, d;
35206cf5b03SBarry Smith   const PetscScalar *xx, *FF;
353c4762a1bSJed Brown   PetscInt           i, M, xs, xm;
354c4762a1bSJed Brown   Vec                xlocal;
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   PetscFunctionBeginUser;
3579566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &xlocal));
358c4762a1bSJed Brown   /*
359c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
360c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
361c4762a1bSJed Brown      By placing code between these two statements, computations can
362c4762a1bSJed Brown      be done while messages are in transition.
363c4762a1bSJed Brown   */
3649566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, x, INSERT_VALUES, xlocal));
3659566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, x, INSERT_VALUES, xlocal));
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown      Get pointers to vector data.
369c4762a1bSJed Brown        - The vector xlocal includes ghost point; the vectors x and f do
370c4762a1bSJed Brown          NOT include ghost points.
371c4762a1bSJed Brown        - Using DMDAVecGetArray() allows accessing the values using global ordering
372c4762a1bSJed Brown   */
37395b2e421SBarry Smith   PetscCall(DMDAVecGetArrayRead(da, xlocal, (void *)&xx));
3749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, f, &ff));
37595b2e421SBarry Smith   PetscCall(DMDAVecGetArrayRead(da, user->F, (void *)&FF));
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   /*
378c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
379c4762a1bSJed Brown        xs, xm  - starting grid index, width of local grid (no ghost points)
380c4762a1bSJed Brown   */
3819566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
3829566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   /*
385c4762a1bSJed Brown      Set function values for boundary points; define local interior grid point range:
386c4762a1bSJed Brown         xsi - starting interior grid index
387c4762a1bSJed Brown         xei - ending interior grid index
388c4762a1bSJed Brown   */
389c4762a1bSJed Brown   if (xs == 0) { /* left boundary */
390c4762a1bSJed Brown     ff[0] = xx[0];
39111486bccSBarry Smith     xs++;
39211486bccSBarry Smith     xm--;
393c4762a1bSJed Brown   }
394c4762a1bSJed Brown   if (xs + xm == M) { /* right boundary */
395c4762a1bSJed Brown     ff[xs + xm - 1] = xx[xs + xm - 1] - 1.0;
396c4762a1bSJed Brown     xm--;
397c4762a1bSJed Brown   }
398c4762a1bSJed Brown 
399c4762a1bSJed Brown   /*
400c4762a1bSJed Brown      Compute function over locally owned part of the grid (interior points only)
401c4762a1bSJed Brown   */
402c4762a1bSJed Brown   d = 1.0 / (user->h * user->h);
403c4762a1bSJed 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];
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   /*
406c4762a1bSJed Brown      Restore vectors
407c4762a1bSJed Brown   */
40895b2e421SBarry Smith   PetscCall(DMDAVecRestoreArrayRead(da, xlocal, (void *)&xx));
4099566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, f, &ff));
41095b2e421SBarry Smith   PetscCall(DMDAVecRestoreArrayRead(da, user->F, (void *)&FF));
4119566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &xlocal));
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413c4762a1bSJed Brown }
414c4762a1bSJed Brown 
415c4762a1bSJed Brown /* ------------------------------------------------------------------- */
416c4762a1bSJed Brown /*
417c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
418c4762a1bSJed Brown 
419c4762a1bSJed Brown    Input Parameters:
420c4762a1bSJed Brown .  snes - the SNES context
421c4762a1bSJed Brown .  x - input vector
422c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
423c4762a1bSJed Brown 
424c4762a1bSJed Brown    Output Parameters:
425c4762a1bSJed Brown .  jac - Jacobian matrix
4267addb90fSBarry Smith .  B - optionally different matrix used to construct the preconditioner
4270b4b7b1cSBarry Smith 
428c4762a1bSJed Brown */
FormJacobian(SNES snes,Vec x,Mat jac,Mat B,PetscCtx ctx)429*2a8381b2SBarry Smith PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, PetscCtx ctx)
430d71ae5a4SJacob Faibussowitsch {
431c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx *)ctx;
432c4762a1bSJed Brown   PetscScalar    *xx, d, A[3];
433c4762a1bSJed Brown   PetscInt        i, j[3], M, xs, xm;
434c4762a1bSJed Brown   DM              da = user->da;
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   PetscFunctionBeginUser;
437c4762a1bSJed Brown   if (user->sjerr) {
4389566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobianDomainError(snes));
4393ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
440c4762a1bSJed Brown   }
441c4762a1bSJed Brown   /*
442c4762a1bSJed Brown      Get pointer to vector data
443c4762a1bSJed Brown   */
4449566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, x, &xx));
4459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
446c4762a1bSJed Brown 
447c4762a1bSJed Brown   /*
448c4762a1bSJed Brown     Get range of locally owned matrix
449c4762a1bSJed Brown   */
4509566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   /*
453c4762a1bSJed Brown      Determine starting and ending local indices for interior grid points.
454c4762a1bSJed Brown      Set Jacobian entries for boundary points.
455c4762a1bSJed Brown   */
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   if (xs == 0) { /* left boundary */
45811486bccSBarry Smith     i    = 0;
45911486bccSBarry Smith     A[0] = 1.0;
460c4762a1bSJed Brown 
4619566063dSJacob Faibussowitsch     PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES));
46211486bccSBarry Smith     xs++;
46311486bccSBarry Smith     xm--;
464c4762a1bSJed Brown   }
465c4762a1bSJed Brown   if (xs + xm == M) { /* right boundary */
466c4762a1bSJed Brown     i    = M - 1;
467c4762a1bSJed Brown     A[0] = 1.0;
4689566063dSJacob Faibussowitsch     PetscCall(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++) {
47911486bccSBarry Smith     j[0] = i - 1;
48011486bccSBarry Smith     j[1] = i;
48111486bccSBarry Smith     j[2] = i + 1;
48211486bccSBarry Smith     A[0] = A[2] = d;
48311486bccSBarry Smith     A[1]        = -2.0 * d + 2.0 * xx[i];
4849566063dSJacob Faibussowitsch     PetscCall(MatSetValues(jac, 1, &i, 3, j, A, INSERT_VALUES));
485c4762a1bSJed Brown   }
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   /*
488c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
489c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
490c4762a1bSJed Brown      By placing code between these two statements, computations can be
491c4762a1bSJed Brown      done while messages are in transition.
492c4762a1bSJed Brown 
493c4762a1bSJed Brown      Also, restore vector.
494c4762a1bSJed Brown   */
495c4762a1bSJed Brown 
4969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
4979566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, x, &xx));
4989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500c4762a1bSJed Brown }
501c4762a1bSJed Brown 
502c4762a1bSJed Brown /* ------------------------------------------------------------------- */
503c4762a1bSJed Brown /*
504c4762a1bSJed Brown    Monitor - Optional user-defined monitoring routine that views the
505c4762a1bSJed Brown    current iterate with an x-window plot. Set by SNESMonitorSet().
506c4762a1bSJed Brown 
507c4762a1bSJed Brown    Input Parameters:
508c4762a1bSJed Brown    snes - the SNES context
509c4762a1bSJed Brown    its - iteration number
510c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
511c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
512c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
513c4762a1bSJed Brown 
514c4762a1bSJed Brown    Note:
515c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
516c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
517c4762a1bSJed Brown  */
Monitor(SNES snes,PetscInt its,PetscReal fnorm,PetscCtx ctx)518*2a8381b2SBarry Smith PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, PetscCtx ctx)
519d71ae5a4SJacob Faibussowitsch {
520c4762a1bSJed Brown   MonitorCtx *monP = (MonitorCtx *)ctx;
521c4762a1bSJed Brown   Vec         x;
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   PetscFunctionBeginUser;
52463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iter = %" PetscInt_FMT ",SNES Function norm %g\n", its, (double)fnorm));
5259566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &x));
5269566063dSJacob Faibussowitsch   PetscCall(VecView(x, monP->viewer));
5273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
528c4762a1bSJed Brown }
529c4762a1bSJed Brown 
530c4762a1bSJed Brown /* ------------------------------------------------------------------- */
531c4762a1bSJed Brown /*
532c4762a1bSJed Brown    PreCheck - Optional user-defined routine that checks the validity of
533c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPreCheck().
534c4762a1bSJed Brown 
535c4762a1bSJed Brown    Input Parameters:
536c4762a1bSJed Brown    snes - the SNES context
537c4762a1bSJed Brown    xcurrent - current solution
538c4762a1bSJed Brown    y - search direction and length
539c4762a1bSJed Brown 
540c4762a1bSJed Brown    Output Parameters:
541c4762a1bSJed Brown    y         - proposed step (search direction and length) (possibly changed)
542c4762a1bSJed Brown    changed_y - tells if the step has changed or not
543c4762a1bSJed Brown  */
PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,PetscBool * changed_y,PetscCtx ctx)544*2a8381b2SBarry Smith PetscErrorCode PreCheck(SNESLineSearch linesearch, Vec xcurrent, Vec y, PetscBool *changed_y, PetscCtx ctx)
545d71ae5a4SJacob Faibussowitsch {
546c4762a1bSJed Brown   PetscFunctionBeginUser;
547c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549c4762a1bSJed Brown }
550c4762a1bSJed Brown 
551c4762a1bSJed Brown /* ------------------------------------------------------------------- */
552c4762a1bSJed Brown /*
553c4762a1bSJed Brown    PostCheck - Optional user-defined routine that checks the validity of
554c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPostCheck().
555c4762a1bSJed Brown 
556c4762a1bSJed Brown    Input Parameters:
557c4762a1bSJed Brown    snes - the SNES context
558c4762a1bSJed Brown    ctx  - optional user-defined context for private data for the
559c4762a1bSJed Brown           monitor routine, as set by SNESLineSearchSetPostCheck()
560c4762a1bSJed Brown    xcurrent - current solution
561c4762a1bSJed Brown    y - search direction and length
562c4762a1bSJed Brown    x    - the new candidate iterate
563c4762a1bSJed Brown 
564c4762a1bSJed Brown    Output Parameters:
565c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
566c4762a1bSJed Brown    x    - current iterate (possibly modified)
567c4762a1bSJed Brown 
568c4762a1bSJed Brown  */
PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool * changed_y,PetscBool * changed_x,PetscCtx ctx)569*2a8381b2SBarry Smith PetscErrorCode PostCheck(SNESLineSearch linesearch, Vec xcurrent, Vec y, Vec x, PetscBool *changed_y, PetscBool *changed_x, PetscCtx ctx)
570d71ae5a4SJacob Faibussowitsch {
571c4762a1bSJed Brown   PetscInt        i, iter, xs, xm;
572c4762a1bSJed Brown   StepCheckCtx   *check;
573c4762a1bSJed Brown   ApplicationCtx *user;
574c4762a1bSJed Brown   PetscScalar    *xa, *xa_last, tmp;
575c4762a1bSJed Brown   PetscReal       rdiff;
576c4762a1bSJed Brown   DM              da;
577c4762a1bSJed Brown   SNES            snes;
578c4762a1bSJed Brown 
579c4762a1bSJed Brown   PetscFunctionBeginUser;
580c4762a1bSJed Brown   *changed_x = PETSC_FALSE;
581c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
582c4762a1bSJed Brown 
5839566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
584c4762a1bSJed Brown   check = (StepCheckCtx *)ctx;
585c4762a1bSJed Brown   user  = check->user;
5869566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   /* iteration 1 indicates we are working on the second iteration */
589c4762a1bSJed Brown   if (iter > 0) {
590c4762a1bSJed Brown     da = user->da;
59163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Checking candidate step at iteration %" PetscInt_FMT " with tolerance %g\n", iter, (double)check->tolerance));
592c4762a1bSJed Brown 
593c4762a1bSJed Brown     /* Access local array data */
5949566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, check->last_step, &xa_last));
5959566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, x, &xa));
5969566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
597c4762a1bSJed Brown 
598c4762a1bSJed Brown     /*
599c4762a1bSJed Brown        If we fail the user-defined check for validity of the candidate iterate,
600c4762a1bSJed Brown        then modify the iterate as we like.  (Note that the particular modification
601c4762a1bSJed Brown        below is intended simply to demonstrate how to manipulate this data, not
602c4762a1bSJed Brown        as a meaningful or appropriate choice.)
603c4762a1bSJed Brown     */
604c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
605c4762a1bSJed Brown       if (!PetscAbsScalar(xa[i])) rdiff = 2 * check->tolerance;
606c4762a1bSJed Brown       else rdiff = PetscAbsScalar((xa[i] - xa_last[i]) / xa[i]);
607c4762a1bSJed Brown       if (rdiff > check->tolerance) {
608c4762a1bSJed Brown         tmp        = xa[i];
609c4762a1bSJed Brown         xa[i]      = .5 * (xa[i] + xa_last[i]);
610c4762a1bSJed Brown         *changed_x = PETSC_TRUE;
61163a3b9bcSJacob 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])));
612c4762a1bSJed Brown       }
613c4762a1bSJed Brown     }
6149566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, check->last_step, &xa_last));
6159566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, x, &xa));
616c4762a1bSJed Brown   }
6179566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, check->last_step));
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
619c4762a1bSJed Brown }
620c4762a1bSJed Brown 
621c4762a1bSJed Brown /* ------------------------------------------------------------------- */
622c4762a1bSJed Brown /*
623c4762a1bSJed Brown    PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used
624c4762a1bSJed Brown    e.g,
625c4762a1bSJed 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
626c4762a1bSJed Brown    Set by SNESLineSearchSetPostCheck().
627c4762a1bSJed Brown 
628c4762a1bSJed Brown    Input Parameters:
629c4762a1bSJed Brown    linesearch - the LineSearch context
630c4762a1bSJed Brown    xcurrent - current solution
631c4762a1bSJed Brown    y - search direction and length
632c4762a1bSJed Brown    x    - the new candidate iterate
633c4762a1bSJed Brown 
634c4762a1bSJed Brown    Output Parameters:
635c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
636c4762a1bSJed Brown    x    - current iterate (possibly modified)
637c4762a1bSJed Brown 
638c4762a1bSJed Brown  */
PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool * changed_y,PetscBool * changed_x,PetscCtx ctx)639*2a8381b2SBarry Smith PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch, Vec xcurrent, Vec y, Vec x, PetscBool *changed_y, PetscBool *changed_x, PetscCtx ctx)
640d71ae5a4SJacob Faibussowitsch {
641c4762a1bSJed Brown   SetSubKSPCtx *check;
642c4762a1bSJed Brown   PetscInt      iter, its, sub_its, maxit;
643c4762a1bSJed Brown   KSP           ksp, sub_ksp, *sub_ksps;
644c4762a1bSJed Brown   PC            pc;
645c4762a1bSJed Brown   PetscReal     ksp_ratio;
646c4762a1bSJed Brown   SNES          snes;
647c4762a1bSJed Brown 
648c4762a1bSJed Brown   PetscFunctionBeginUser;
6499566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
650c4762a1bSJed Brown   check = (SetSubKSPCtx *)ctx;
6519566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
6529566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
6539566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
6549566063dSJacob Faibussowitsch   PetscCall(PCBJacobiGetSubKSP(pc, NULL, NULL, &sub_ksps));
655c4762a1bSJed Brown   sub_ksp = sub_ksps[0];
6569566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(ksp, &its));         /* outer KSP iteration number */
6579566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(sub_ksp, &sub_its)); /* inner KSP iteration number */
658c4762a1bSJed Brown 
659c4762a1bSJed Brown   if (iter) {
66063a3b9bcSJacob 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));
66157508eceSPierre Jolivet     ksp_ratio = (PetscReal)its / check->its0;
662c4762a1bSJed Brown     maxit     = (PetscInt)(ksp_ratio * sub_its + 0.5);
663c4762a1bSJed Brown     if (maxit < 2) maxit = 2;
664fb842aefSJose E. Roman     PetscCall(KSPSetTolerances(sub_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, maxit));
66563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "    ...ksp_ratio %g, new maxit %" PetscInt_FMT "\n\n", (double)ksp_ratio, maxit));
666c4762a1bSJed Brown   }
667c4762a1bSJed Brown   check->its0 = its; /* save current outer KSP iteration number */
6683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
669c4762a1bSJed Brown }
670c4762a1bSJed Brown 
671c4762a1bSJed Brown /* ------------------------------------------------------------------- */
672c4762a1bSJed Brown /*
673c4762a1bSJed Brown    MatrixFreePreconditioner - This routine demonstrates the use of a
674c4762a1bSJed Brown    user-provided preconditioner.  This code implements just the null
675c4762a1bSJed Brown    preconditioner, which of course is not recommended for general use.
676c4762a1bSJed Brown 
677c4762a1bSJed Brown    Input Parameters:
678c4762a1bSJed Brown +  pc - preconditioner
679c4762a1bSJed Brown -  x - input vector
680c4762a1bSJed Brown 
681c4762a1bSJed Brown    Output Parameter:
682c4762a1bSJed Brown .  y - preconditioned vector
683c4762a1bSJed Brown */
MatrixFreePreconditioner(PC pc,Vec x,Vec y)684d71ae5a4SJacob Faibussowitsch PetscErrorCode MatrixFreePreconditioner(PC pc, Vec x, Vec y)
685d71ae5a4SJacob Faibussowitsch {
6863ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
6879566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, y));
6883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
689c4762a1bSJed Brown }
690c4762a1bSJed Brown 
691c4762a1bSJed Brown /*TEST
692c4762a1bSJed Brown 
693c4762a1bSJed Brown    test:
694c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
695c4762a1bSJed Brown 
696c4762a1bSJed Brown    test:
697c4762a1bSJed Brown       suffix: 2
698c4762a1bSJed Brown       nsize: 3
699c4762a1bSJed Brown       args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
700c4762a1bSJed Brown 
701c4762a1bSJed Brown    test:
702c4762a1bSJed Brown       suffix: 3
703c4762a1bSJed Brown       nsize: 2
704c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
705c4762a1bSJed Brown 
706c4762a1bSJed Brown    test:
707c4762a1bSJed Brown       suffix: 4
708c4762a1bSJed Brown       args: -nox -pre_check_iterates -post_check_iterates
709c4762a1bSJed Brown 
710c4762a1bSJed Brown    test:
711c4762a1bSJed Brown       suffix: 5
712c4762a1bSJed Brown       requires: double !complex !single
713c4762a1bSJed Brown       nsize: 2
714c4762a1bSJed Brown       args: -nox -snes_test_jacobian -snes_test_jacobian_view
715c4762a1bSJed Brown 
716c4762a1bSJed Brown    test:
717c4762a1bSJed Brown       suffix: 6
718c4762a1bSJed Brown       requires: double !complex !single
719c4762a1bSJed Brown       nsize: 4
720c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1
721c4762a1bSJed Brown 
722c4762a1bSJed Brown    test:
723c4762a1bSJed Brown       suffix: 7
724c4762a1bSJed Brown       requires: double !complex !single
725c4762a1bSJed Brown       nsize: 4
726c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1
727c4762a1bSJed Brown 
728c4762a1bSJed Brown    test:
729c4762a1bSJed Brown       suffix: 8
730c4762a1bSJed Brown       requires: double !complex !single
731c4762a1bSJed Brown       nsize: 4
732c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1
733c4762a1bSJed Brown 
734c4762a1bSJed Brown    test:
735c4762a1bSJed Brown       suffix: 9
736c4762a1bSJed Brown       requires: double !complex !single
737c4762a1bSJed Brown       nsize: 4
738c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1
739c4762a1bSJed Brown 
740c4762a1bSJed Brown    test:
741c4762a1bSJed Brown       suffix: 10
742c4762a1bSJed Brown       requires: double !complex !single
743c4762a1bSJed Brown       nsize: 4
744c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1
745c4762a1bSJed Brown 
746c4762a1bSJed Brown    test:
747c4762a1bSJed Brown       suffix: 11
748c4762a1bSJed Brown       requires: double !complex !single
749c4762a1bSJed Brown       nsize: 4
75057715debSLisandro 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
751c4762a1bSJed Brown 
752c0558f20SBarry Smith    test:
753c0558f20SBarry Smith       suffix: 12
754c0558f20SBarry Smith       args: -view_initial
755c0558f20SBarry Smith       filter: grep -v "type:"
756c0558f20SBarry Smith 
75741ba4c6cSHeeho Park    test:
75841ba4c6cSHeeho Park       suffix: 13
75941ba4c6cSHeeho Park       requires: double !complex !single
76041ba4c6cSHeeho Park       nsize: 4
76141ba4c6cSHeeho Park       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1
76241ba4c6cSHeeho Park 
763c4762a1bSJed Brown TEST*/
764