xref: /petsc/src/snes/tutorials/ex2.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
2c4762a1bSJed Brown This example employs a user-defined monitoring routine.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Include "petscdraw.h" so that we can use PETSc drawing routines.
6c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
7c4762a1bSJed Brown    file automatically includes:
8c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
9c4762a1bSJed Brown      petscmat.h - matrices
10c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
11c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
12c4762a1bSJed Brown      petscksp.h   - linear solvers
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petscsnes.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /*
18c4762a1bSJed Brown    User-defined routines
19c4762a1bSJed Brown */
20c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
21c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
22c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
23c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *);
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    User-defined context for monitoring
27c4762a1bSJed Brown */
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown   PetscViewer viewer;
30c4762a1bSJed Brown } MonitorCtx;
31c4762a1bSJed Brown 
main(int argc,char ** argv)32d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
33d71ae5a4SJacob Faibussowitsch {
34c4762a1bSJed Brown   SNES        snes;       /* SNES context */
35c4762a1bSJed Brown   Vec         x, r, F, U; /* vectors */
36c4762a1bSJed Brown   Mat         J;          /* Jacobian matrix */
37c4762a1bSJed Brown   MonitorCtx  monP;       /* monitoring context */
38c4762a1bSJed Brown   PetscInt    its, n = 5, i, maxit, maxf;
39c4762a1bSJed Brown   PetscMPIInt size;
40c4762a1bSJed Brown   PetscScalar h, xp, v, none = -1.0;
41c4762a1bSJed Brown   PetscReal   abstol, rtol, stol, norm;
42c4762a1bSJed Brown 
43327415f7SBarry Smith   PetscFunctionBeginUser;
44c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
46be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
48c4762a1bSJed Brown   h = 1.0 / (n - 1);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Create nonlinear solver context
52c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
58c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59c4762a1bSJed Brown   /*
60c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
61c4762a1bSJed Brown   */
629566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
639566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, n));
649566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &F));
679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /*
70c4762a1bSJed Brown      Set function evaluation routine and vector
71c4762a1bSJed Brown   */
729566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
76c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, n, n));
809566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
819566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /*
84c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
85c4762a1bSJed Brown      routine. User can override with:
86c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
87c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
88c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
897addb90fSBarry Smith      -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
90c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
91c4762a1bSJed Brown                          products within Newton-Krylov method
92c4762a1bSJed Brown   */
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
98c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown      Set an optional user-defined monitoring routine
102c4762a1bSJed Brown   */
1039566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, 0, 0, 0, 400, 400, &monP.viewer));
1049566063dSJacob Faibussowitsch   PetscCall(SNESMonitorSet(snes, Monitor, &monP, 0));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
108c4762a1bSJed Brown   */
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution"));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution"));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*
113c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
114c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
115c4762a1bSJed Brown   */
1169566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /*
119c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
120c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
121c4762a1bSJed Brown      the option -snes_view
122c4762a1bSJed Brown   */
1239566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf));
12463a3b9bcSJacob 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));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Initialize application:
128dd8e379bSPierre Jolivet      Store right-hand side of PDE and exact solution
129c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   xp = 0.0;
132c4762a1bSJed Brown   for (i = 0; i < n; i++) {
133c4762a1bSJed Brown     v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
1349566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES));
135c4762a1bSJed Brown     v = xp * xp * xp;
1369566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES));
137c4762a1bSJed Brown     xp += h;
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
142c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143c4762a1bSJed Brown   /*
144c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
145c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
146c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
147c4762a1bSJed Brown      this vector to zero by calling VecSet().
148c4762a1bSJed Brown   */
1499566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1509566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1519566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
15263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155c4762a1bSJed Brown      Check solution and clean up
156c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /*
159c4762a1bSJed Brown      Check the error
160c4762a1bSJed Brown   */
1619566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, none, U));
1629566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
16363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, Iterations %" PetscInt_FMT "\n", (double)norm, its));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /*
166c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
167c4762a1bSJed Brown      are no longer needed.
168c4762a1bSJed Brown   */
1699371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1709371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1719371c9d4SSatish Balay   PetscCall(VecDestroy(&U));
1729371c9d4SSatish Balay   PetscCall(VecDestroy(&F));
1739371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1749371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1759566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&monP.viewer));
1769566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
177b122ec5aSJacob Faibussowitsch   return 0;
178c4762a1bSJed Brown }
179c4762a1bSJed Brown /* ------------------------------------------------------------------- */
180c4762a1bSJed Brown /*
181c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
182c4762a1bSJed Brown 
183c4762a1bSJed Brown    Input/Output Parameter:
184c4762a1bSJed Brown .  x - the solution vector
185c4762a1bSJed Brown */
FormInitialGuess(Vec x)186d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
187d71ae5a4SJacob Faibussowitsch {
1883ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1894d86920dSPierre Jolivet   PetscCall(VecSet(x, 0.5));
1903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown /* ------------------------------------------------------------------- */
193c4762a1bSJed Brown /*
194c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
195c4762a1bSJed Brown 
196c4762a1bSJed Brown    Input Parameters:
197c4762a1bSJed Brown .  snes - the SNES context
198c4762a1bSJed Brown .  x - input vector
199c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
200c4762a1bSJed Brown 
201c4762a1bSJed Brown    Output Parameter:
202c4762a1bSJed Brown .  f - function vector
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    Note:
205c4762a1bSJed Brown    The user-defined context can contain any application-specific data
206c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
207c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
208dd8e379bSPierre Jolivet    a vector containing the right-hand side of the discretized PDE.
209c4762a1bSJed Brown  */
210c4762a1bSJed Brown 
FormFunction(SNES snes,Vec x,Vec f,PetscCtx ctx)211*2a8381b2SBarry Smith PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, PetscCtx ctx)
212d71ae5a4SJacob Faibussowitsch {
213c4762a1bSJed Brown   Vec                g = (Vec)ctx;
214c4762a1bSJed Brown   const PetscScalar *xx, *gg;
215c4762a1bSJed Brown   PetscScalar       *ff, d;
216c4762a1bSJed Brown   PetscInt           i, n;
217c4762a1bSJed Brown 
2183ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
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));
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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
2597addb90fSBarry Smith .  B - optionally different matrix used to construct the preconditioner
260c4762a1bSJed Brown 
261c4762a1bSJed Brown */
262c4762a1bSJed Brown 
FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void * dummy)263d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
264d71ae5a4SJacob Faibussowitsch {
265c4762a1bSJed Brown   const PetscScalar *xx;
266c4762a1bSJed Brown   PetscScalar        A[3], d;
267c4762a1bSJed Brown   PetscInt           i, n, j[3];
268c4762a1bSJed Brown 
2693ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
270c4762a1bSJed Brown   /*
271c4762a1bSJed Brown      Get pointer to vector data
272c4762a1bSJed Brown   */
2739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /*
276c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
277c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
278c4762a1bSJed Brown         row at once.
279c4762a1bSJed Brown   */
2809566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2819371c9d4SSatish Balay   d = (PetscReal)(n - 1);
2829371c9d4SSatish Balay   d = d * d;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /*
285c4762a1bSJed Brown      Interior grid points
286c4762a1bSJed Brown   */
287c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
2889371c9d4SSatish Balay     j[0] = i - 1;
2899371c9d4SSatish Balay     j[1] = i;
2909371c9d4SSatish Balay     j[2] = i + 1;
2919371c9d4SSatish Balay     A[0] = A[2] = d;
2929371c9d4SSatish Balay     A[1]        = -2.0 * d + 2.0 * xx[i];
2939566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
294c4762a1bSJed Brown   }
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   /*
297c4762a1bSJed Brown      Boundary points
298c4762a1bSJed Brown   */
2999371c9d4SSatish Balay   i    = 0;
3009371c9d4SSatish Balay   A[0] = 1.0;
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
303c4762a1bSJed Brown 
3049371c9d4SSatish Balay   i    = n - 1;
3059371c9d4SSatish Balay   A[0] = 1.0;
306c4762a1bSJed Brown 
3079566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   /*
310c4762a1bSJed Brown      Restore vector
311c4762a1bSJed Brown   */
3129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /*
315c4762a1bSJed Brown      Assemble matrix
316c4762a1bSJed Brown   */
3179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
319c4762a1bSJed Brown   if (jac != B) {
3209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
3219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
322c4762a1bSJed Brown   }
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
324c4762a1bSJed Brown }
325c4762a1bSJed Brown /* ------------------------------------------------------------------- */
326c4762a1bSJed Brown /*
327c4762a1bSJed Brown    Monitor - User-defined monitoring routine that views the
328c4762a1bSJed Brown    current iterate with an x-window plot.
329c4762a1bSJed Brown 
330c4762a1bSJed Brown    Input Parameters:
331c4762a1bSJed Brown    snes - the SNES context
332c4762a1bSJed Brown    its - iteration number
333c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
334c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
335c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    Note:
338c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
339c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
340c4762a1bSJed Brown  */
Monitor(SNES snes,PetscInt its,PetscReal fnorm,PetscCtx ctx)341*2a8381b2SBarry Smith PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, PetscCtx ctx)
342d71ae5a4SJacob Faibussowitsch {
343c4762a1bSJed Brown   MonitorCtx         *monP = (MonitorCtx *)ctx;
344c4762a1bSJed Brown   Vec                 x;
3452d157150SStefano Zampini   SNESConvergedReason reason;
346c4762a1bSJed Brown 
3473ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
34863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iter = %" PetscInt_FMT ", SNES Function norm %g\n", its, (double)fnorm));
3492d157150SStefano Zampini   PetscCall(SNESGetConvergedReason(snes, &reason));
3509566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &x));
3519566063dSJacob Faibussowitsch   PetscCall(VecView(x, monP->viewer));
3522d157150SStefano Zampini   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  converged = %s\n", SNESConvergedReasons[reason]));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354c4762a1bSJed Brown }
355c4762a1bSJed Brown 
356c4762a1bSJed Brown /*TEST
357c4762a1bSJed Brown 
358c4762a1bSJed Brown    test:
359c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
360c4762a1bSJed Brown 
361c4762a1bSJed Brown    test:
362c4762a1bSJed Brown       suffix: 2
363c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontr -snes_view
364c4762a1bSJed Brown       requires: !single
365c4762a1bSJed Brown 
366c4762a1bSJed Brown    test:
367c4762a1bSJed Brown       suffix: 3
368c4762a1bSJed 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
369c4762a1bSJed Brown 
37041ba4c6cSHeeho Park    test:
37141ba4c6cSHeeho Park       suffix: 4
37241ba4c6cSHeeho Park       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontrdc -snes_view
37341ba4c6cSHeeho Park       requires: !single
37441ba4c6cSHeeho Park 
3752d157150SStefano Zampini    test:
3762d157150SStefano Zampini       suffix: 5
3777a14dab9SStefano Zampini       filter: grep -v atol | sed -e "s/CONVERGED_ITS/DIVERGED_MAX_IT/g" | sed -e "s/CONVERGED_FNORM_RELATIVE/DIVERGED_MAX_IT/g"
3787a14dab9SStefano Zampini       args: -nox -snes_type {{newtonls newtontr ncg ngmres qn anderson nrichardson ms ksponly ksptransposeonly vinewtonrsls vinewtonssls fas ms}} -snes_max_it 1
3792d157150SStefano Zampini       requires: !single
3802d157150SStefano Zampini 
381c4762a1bSJed Brown TEST*/
382