xref: /petsc/src/snes/tutorials/ex2.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3c4762a1bSJed Brown This example employs a user-defined monitoring routine.\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Include "petscdraw.h" so that we can use PETSc drawing routines.
7c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
8c4762a1bSJed Brown    file automatically includes:
9c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
10c4762a1bSJed Brown      petscmat.h - matrices
11c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
12c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
13c4762a1bSJed Brown      petscksp.h   - linear solvers
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown #include <petscsnes.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    User-defined routines
20c4762a1bSJed Brown */
21c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
22c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
23c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
24c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown    User-defined context for monitoring
28c4762a1bSJed Brown */
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   PetscViewer viewer;
31c4762a1bSJed Brown } MonitorCtx;
32c4762a1bSJed Brown 
33*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
34*d71ae5a4SJacob Faibussowitsch {
35c4762a1bSJed Brown   SNES        snes;       /* SNES context */
36c4762a1bSJed Brown   Vec         x, r, F, U; /* vectors */
37c4762a1bSJed Brown   Mat         J;          /* Jacobian matrix */
38c4762a1bSJed Brown   MonitorCtx  monP;       /* monitoring context */
39c4762a1bSJed Brown   PetscInt    its, n = 5, i, maxit, maxf;
40c4762a1bSJed Brown   PetscMPIInt size;
41c4762a1bSJed Brown   PetscScalar h, xp, v, none = -1.0;
42c4762a1bSJed Brown   PetscReal   abstol, rtol, stol, norm;
43c4762a1bSJed Brown 
44327415f7SBarry Smith   PetscFunctionBeginUser;
459566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
47be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
49c4762a1bSJed Brown   h = 1.0 / (n - 1);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52c4762a1bSJed Brown      Create nonlinear solver context
53c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60c4762a1bSJed Brown   /*
61c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
62c4762a1bSJed Brown   */
639566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
649566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, n));
659566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &F));
689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /*
71c4762a1bSJed Brown      Set function evaluation routine and vector
72c4762a1bSJed Brown   */
739566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
77c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
809566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, n, n));
819566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
829566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /*
85c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
86c4762a1bSJed Brown      routine. User can override with:
87c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
88c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
89c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
90c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
91c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
92c4762a1bSJed Brown                          products within Newton-Krylov method
93c4762a1bSJed Brown   */
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
99c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      Set an optional user-defined monitoring routine
103c4762a1bSJed Brown   */
1049566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, 0, 0, 0, 400, 400, &monP.viewer));
1059566063dSJacob Faibussowitsch   PetscCall(SNESMonitorSet(snes, Monitor, &monP, 0));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /*
108c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
109c4762a1bSJed Brown   */
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution"));
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution"));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /*
114c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
115c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
116c4762a1bSJed Brown   */
1179566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /*
120c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
121c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
122c4762a1bSJed Brown      the option -snes_view
123c4762a1bSJed Brown   */
1249566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf));
12563a3b9bcSJacob 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));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128c4762a1bSJed Brown      Initialize application:
129c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
130c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   xp = 0.0;
133c4762a1bSJed Brown   for (i = 0; i < n; i++) {
134c4762a1bSJed Brown     v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
1359566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES));
136c4762a1bSJed Brown     v = xp * xp * xp;
1379566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES));
138c4762a1bSJed Brown     xp += h;
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
143c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144c4762a1bSJed Brown   /*
145c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
146c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
147c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
148c4762a1bSJed Brown      this vector to zero by calling VecSet().
149c4762a1bSJed Brown   */
1509566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1519566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1529566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
15363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156c4762a1bSJed Brown      Check solution and clean up
157c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /*
160c4762a1bSJed Brown      Check the error
161c4762a1bSJed Brown   */
1629566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, none, U));
1639566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
16463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, Iterations %" PetscInt_FMT "\n", (double)norm, its));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /*
167c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
168c4762a1bSJed Brown      are no longer needed.
169c4762a1bSJed Brown   */
1709371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1719371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1729371c9d4SSatish Balay   PetscCall(VecDestroy(&U));
1739371c9d4SSatish Balay   PetscCall(VecDestroy(&F));
1749371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1759371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1769566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&monP.viewer));
1779566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
178b122ec5aSJacob Faibussowitsch   return 0;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown /* ------------------------------------------------------------------- */
181c4762a1bSJed Brown /*
182c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
183c4762a1bSJed Brown 
184c4762a1bSJed Brown    Input/Output Parameter:
185c4762a1bSJed Brown .  x - the solution vector
186c4762a1bSJed Brown */
187*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
188*d71ae5a4SJacob Faibussowitsch {
189c4762a1bSJed Brown   PetscScalar pfive = .50;
1909566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
191c4762a1bSJed Brown   return 0;
192c4762a1bSJed Brown }
193c4762a1bSJed Brown /* ------------------------------------------------------------------- */
194c4762a1bSJed Brown /*
195c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
196c4762a1bSJed Brown 
197c4762a1bSJed Brown    Input Parameters:
198c4762a1bSJed Brown .  snes - the SNES context
199c4762a1bSJed Brown .  x - input vector
200c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
201c4762a1bSJed Brown 
202c4762a1bSJed Brown    Output Parameter:
203c4762a1bSJed Brown .  f - function vector
204c4762a1bSJed Brown 
205c4762a1bSJed Brown    Note:
206c4762a1bSJed Brown    The user-defined context can contain any application-specific data
207c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
208c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
209c4762a1bSJed Brown    a vector containing the right-hand-side of the discretized PDE.
210c4762a1bSJed Brown  */
211c4762a1bSJed Brown 
212*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx)
213*d71ae5a4SJacob Faibussowitsch {
214c4762a1bSJed Brown   Vec                g = (Vec)ctx;
215c4762a1bSJed Brown   const PetscScalar *xx, *gg;
216c4762a1bSJed Brown   PetscScalar       *ff, d;
217c4762a1bSJed Brown   PetscInt           i, n;
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /*
220c4762a1bSJed Brown      Get pointers to vector data.
221c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
222c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
223c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
224c4762a1bSJed Brown          the array.
225c4762a1bSJed Brown   */
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
2289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g, &gg));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /*
231c4762a1bSJed Brown      Compute function
232c4762a1bSJed Brown   */
2339566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2349371c9d4SSatish Balay   d     = (PetscReal)(n - 1);
2359371c9d4SSatish Balay   d     = d * d;
236c4762a1bSJed Brown   ff[0] = xx[0];
237c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) ff[i] = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - gg[i];
238c4762a1bSJed Brown   ff[n - 1] = xx[n - 1] - 1.0;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   /*
241c4762a1bSJed Brown      Restore vectors
242c4762a1bSJed Brown   */
2439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
2459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g, &gg));
246c4762a1bSJed Brown   return 0;
247c4762a1bSJed Brown }
248c4762a1bSJed Brown /* ------------------------------------------------------------------- */
249c4762a1bSJed Brown /*
250c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    Input Parameters:
253c4762a1bSJed Brown .  snes - the SNES context
254c4762a1bSJed Brown .  x - input vector
255c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
256c4762a1bSJed Brown 
257c4762a1bSJed Brown    Output Parameters:
258c4762a1bSJed Brown .  jac - Jacobian matrix
259c4762a1bSJed Brown .  B - optionally different preconditioning matrix
260c4762a1bSJed Brown 
261c4762a1bSJed Brown */
262c4762a1bSJed Brown 
263*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
264*d71ae5a4SJacob Faibussowitsch {
265c4762a1bSJed Brown   const PetscScalar *xx;
266c4762a1bSJed Brown   PetscScalar        A[3], d;
267c4762a1bSJed Brown   PetscInt           i, n, j[3];
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /*
270c4762a1bSJed Brown      Get pointer to vector data
271c4762a1bSJed Brown   */
2729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /*
275c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
276c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
277c4762a1bSJed Brown         row at once.
278c4762a1bSJed Brown   */
2799566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2809371c9d4SSatish Balay   d = (PetscReal)(n - 1);
2819371c9d4SSatish Balay   d = d * d;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /*
284c4762a1bSJed Brown      Interior grid points
285c4762a1bSJed Brown   */
286c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
2879371c9d4SSatish Balay     j[0] = i - 1;
2889371c9d4SSatish Balay     j[1] = i;
2899371c9d4SSatish Balay     j[2] = i + 1;
2909371c9d4SSatish Balay     A[0] = A[2] = d;
2919371c9d4SSatish Balay     A[1]        = -2.0 * d + 2.0 * xx[i];
2929566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /*
296c4762a1bSJed Brown      Boundary points
297c4762a1bSJed Brown   */
2989371c9d4SSatish Balay   i    = 0;
2999371c9d4SSatish Balay   A[0] = 1.0;
300c4762a1bSJed Brown 
3019566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
302c4762a1bSJed Brown 
3039371c9d4SSatish Balay   i    = n - 1;
3049371c9d4SSatish Balay   A[0] = 1.0;
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /*
309c4762a1bSJed Brown      Restore vector
310c4762a1bSJed Brown   */
3119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /*
314c4762a1bSJed Brown      Assemble matrix
315c4762a1bSJed Brown   */
3169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
318c4762a1bSJed Brown   if (jac != B) {
3199566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
3209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
321c4762a1bSJed Brown   }
322c4762a1bSJed Brown   return 0;
323c4762a1bSJed Brown }
324c4762a1bSJed Brown /* ------------------------------------------------------------------- */
325c4762a1bSJed Brown /*
326c4762a1bSJed Brown    Monitor - User-defined monitoring routine that views the
327c4762a1bSJed Brown    current iterate with an x-window plot.
328c4762a1bSJed Brown 
329c4762a1bSJed Brown    Input Parameters:
330c4762a1bSJed Brown    snes - the SNES context
331c4762a1bSJed Brown    its - iteration number
332c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
333c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
334c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    Note:
337c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
338c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
339c4762a1bSJed Brown  */
340*d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, void *ctx)
341*d71ae5a4SJacob Faibussowitsch {
342c4762a1bSJed Brown   MonitorCtx *monP = (MonitorCtx *)ctx;
343c4762a1bSJed Brown   Vec         x;
344c4762a1bSJed Brown 
34563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iter = %" PetscInt_FMT ", SNES Function norm %g\n", its, (double)fnorm));
3469566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &x));
3479566063dSJacob Faibussowitsch   PetscCall(VecView(x, monP->viewer));
348c4762a1bSJed Brown   return 0;
349c4762a1bSJed Brown }
350c4762a1bSJed Brown 
351c4762a1bSJed Brown /*TEST
352c4762a1bSJed Brown 
353c4762a1bSJed Brown    test:
354c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
355c4762a1bSJed Brown 
356c4762a1bSJed Brown    test:
357c4762a1bSJed Brown       suffix: 2
358c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontr -snes_view
359c4762a1bSJed Brown       requires: !single
360c4762a1bSJed Brown 
361c4762a1bSJed Brown    test:
362c4762a1bSJed Brown       suffix: 3
363c4762a1bSJed Brown       args: -nox -malloc no -options_left no -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
364c4762a1bSJed Brown 
36541ba4c6cSHeeho Park    test:
36641ba4c6cSHeeho Park       suffix: 4
36741ba4c6cSHeeho Park       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontrdc -snes_view
36841ba4c6cSHeeho Park       requires: !single
36941ba4c6cSHeeho Park 
370c4762a1bSJed Brown TEST*/
371