xref: /petsc/src/snes/tutorials/ex2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 /*T
6c4762a1bSJed Brown    Concepts: SNES^basic uniprocessor example
7c4762a1bSJed Brown    Concepts: SNES^setting a user-defined monitoring routine
8c4762a1bSJed Brown    Processors: 1
9c4762a1bSJed Brown T*/
10c4762a1bSJed Brown 
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown    Include "petscdraw.h" so that we can use PETSc drawing routines.
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   petscvec.h - vectors
16c4762a1bSJed Brown      petscmat.h - matrices
17c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
18c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
19c4762a1bSJed Brown      petscksp.h   - linear solvers
20c4762a1bSJed Brown */
21c4762a1bSJed Brown 
22c4762a1bSJed Brown #include <petscsnes.h>
23c4762a1bSJed Brown 
24c4762a1bSJed Brown /*
25c4762a1bSJed Brown    User-defined routines
26c4762a1bSJed Brown */
27c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
28c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
29c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
30c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    User-defined context for monitoring
34c4762a1bSJed Brown */
35c4762a1bSJed Brown typedef struct {
36c4762a1bSJed Brown   PetscViewer viewer;
37c4762a1bSJed Brown } MonitorCtx;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown int main(int argc,char **argv)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   SNES           snes;                   /* SNES context */
42c4762a1bSJed Brown   Vec            x,r,F,U;             /* vectors */
43c4762a1bSJed Brown   Mat            J;                      /* Jacobian matrix */
44c4762a1bSJed Brown   MonitorCtx     monP;                   /* monitoring context */
45c4762a1bSJed Brown   PetscInt       its,n = 5,i,maxit,maxf;
46c4762a1bSJed Brown   PetscMPIInt    size;
47c4762a1bSJed Brown   PetscScalar    h,xp,v,none = -1.0;
48c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
49c4762a1bSJed Brown 
50*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
515f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
522c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!");
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
54c4762a1bSJed Brown   h    = 1.0/(n-1);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Create nonlinear solver context
58c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59c4762a1bSJed Brown 
605f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
64c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65c4762a1bSJed Brown   /*
66c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
67c4762a1bSJed Brown   */
685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,n));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /*
76c4762a1bSJed Brown      Set function evaluation routine and vector
77c4762a1bSJed Brown   */
785f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormFunction,(void*)F));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
82c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83c4762a1bSJed Brown 
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,3,NULL));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /*
90c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
91c4762a1bSJed Brown      routine. User can override with:
92c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
93c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
94c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
95c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
96c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
97c4762a1bSJed Brown                          products within Newton-Krylov method
98c4762a1bSJed Brown   */
99c4762a1bSJed Brown 
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
104c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Set an optional user-defined monitoring routine
108c4762a1bSJed Brown   */
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorSet(snes,Monitor,&monP,0));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*
113c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
114c4762a1bSJed Brown   */
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution"));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /*
119c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
120c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
121c4762a1bSJed Brown   */
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /*
125c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
126c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
127c4762a1bSJed Brown      the option -snes_view
128c4762a1bSJed Brown   */
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown      Initialize application:
134c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
135c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   xp = 0.0;
138c4762a1bSJed Brown   for (i=0; i<n; i++) {
139c4762a1bSJed Brown     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES));
141c4762a1bSJed Brown     v    = xp*xp*xp;
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(U,1,&i,&v,INSERT_VALUES));
143c4762a1bSJed Brown     xp  += h;
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
148c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149c4762a1bSJed Brown   /*
150c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
151c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
152c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
153c4762a1bSJed Brown      this vector to zero by calling VecSet().
154c4762a1bSJed Brown   */
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %D\n\n",its));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161c4762a1bSJed Brown      Check solution and clean up
162c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /*
165c4762a1bSJed Brown      Check the error
166c4762a1bSJed Brown   */
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,none,U));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /*
172c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
173c4762a1bSJed Brown      are no longer needed.
174c4762a1bSJed Brown   */
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));  CHKERRQ(VecDestroy(&r));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));  CHKERRQ(VecDestroy(&F));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));  CHKERRQ(SNESDestroy(&snes));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&monP.viewer));
179*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
180*b122ec5aSJacob Faibussowitsch   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown /* ------------------------------------------------------------------- */
183c4762a1bSJed Brown /*
184c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Input/Output Parameter:
187c4762a1bSJed Brown .  x - the solution vector
188c4762a1bSJed Brown */
189c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
190c4762a1bSJed Brown {
191c4762a1bSJed Brown   PetscScalar    pfive = .50;
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
193c4762a1bSJed Brown   return 0;
194c4762a1bSJed Brown }
195c4762a1bSJed Brown /* ------------------------------------------------------------------- */
196c4762a1bSJed Brown /*
197c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    Input Parameters:
200c4762a1bSJed Brown .  snes - the SNES context
201c4762a1bSJed Brown .  x - input vector
202c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    Output Parameter:
205c4762a1bSJed Brown .  f - function vector
206c4762a1bSJed Brown 
207c4762a1bSJed Brown    Note:
208c4762a1bSJed Brown    The user-defined context can contain any application-specific data
209c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
210c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
211c4762a1bSJed Brown    a vector containing the right-hand-side of the discretized PDE.
212c4762a1bSJed Brown  */
213c4762a1bSJed Brown 
214c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   Vec               g = (Vec)ctx;
217c4762a1bSJed Brown   const PetscScalar *xx,*gg;
218c4762a1bSJed Brown   PetscScalar       *ff,d;
219c4762a1bSJed Brown   PetscInt          i,n;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /*
222c4762a1bSJed Brown      Get pointers to vector data.
223c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
224c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
225c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
226c4762a1bSJed Brown          the array.
227c4762a1bSJed Brown   */
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&ff));
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(g,&gg));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown      Compute function
234c4762a1bSJed Brown   */
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
236c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
237c4762a1bSJed Brown   ff[0] = xx[0];
238c4762a1bSJed 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];
239c4762a1bSJed Brown   ff[n-1] = xx[n-1] - 1.0;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /*
242c4762a1bSJed Brown      Restore vectors
243c4762a1bSJed Brown   */
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&ff));
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(g,&gg));
247c4762a1bSJed Brown   return 0;
248c4762a1bSJed Brown }
249c4762a1bSJed Brown /* ------------------------------------------------------------------- */
250c4762a1bSJed Brown /*
251c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
252c4762a1bSJed Brown 
253c4762a1bSJed Brown    Input Parameters:
254c4762a1bSJed Brown .  snes - the SNES context
255c4762a1bSJed Brown .  x - input vector
256c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    Output Parameters:
259c4762a1bSJed Brown .  jac - Jacobian matrix
260c4762a1bSJed Brown .  B - optionally different preconditioning matrix
261c4762a1bSJed Brown 
262c4762a1bSJed Brown */
263c4762a1bSJed Brown 
264c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
265c4762a1bSJed Brown {
266c4762a1bSJed Brown   const PetscScalar *xx;
267c4762a1bSJed Brown   PetscScalar       A[3],d;
268c4762a1bSJed Brown   PetscInt          i,n,j[3];
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   /*
271c4762a1bSJed Brown      Get pointer to vector data
272c4762a1bSJed Brown   */
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(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   */
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&n));
281c4762a1bSJed Brown   d    = (PetscReal)(n - 1); d = d*d;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /*
284c4762a1bSJed Brown      Interior grid points
285c4762a1bSJed Brown   */
286c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
287c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
288c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
2895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
290c4762a1bSJed Brown   }
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /*
293c4762a1bSJed Brown      Boundary points
294c4762a1bSJed Brown   */
295c4762a1bSJed Brown   i = 0;   A[0] = 1.0;
296c4762a1bSJed Brown 
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   i = n-1; A[0] = 1.0;
300c4762a1bSJed Brown 
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /*
304c4762a1bSJed Brown      Restore vector
305c4762a1bSJed Brown   */
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /*
309c4762a1bSJed Brown      Assemble matrix
310c4762a1bSJed Brown   */
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
313c4762a1bSJed Brown   if (jac != B) {
3145f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
3155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
316c4762a1bSJed Brown   }
317c4762a1bSJed Brown   return 0;
318c4762a1bSJed Brown }
319c4762a1bSJed Brown /* ------------------------------------------------------------------- */
320c4762a1bSJed Brown /*
321c4762a1bSJed Brown    Monitor - User-defined monitoring routine that views the
322c4762a1bSJed Brown    current iterate with an x-window plot.
323c4762a1bSJed Brown 
324c4762a1bSJed Brown    Input Parameters:
325c4762a1bSJed Brown    snes - the SNES context
326c4762a1bSJed Brown    its - iteration number
327c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
328c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
329c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
330c4762a1bSJed Brown 
331c4762a1bSJed Brown    Note:
332c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
333c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
334c4762a1bSJed Brown  */
335c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
336c4762a1bSJed Brown {
337c4762a1bSJed Brown   MonitorCtx     *monP = (MonitorCtx*) ctx;
338c4762a1bSJed Brown   Vec            x;
339c4762a1bSJed Brown 
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"iter = %D, SNES Function norm %g\n",its,(double)fnorm));
3415f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes,&x));
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,monP->viewer));
343c4762a1bSJed Brown   return 0;
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
346c4762a1bSJed Brown /*TEST
347c4762a1bSJed Brown 
348c4762a1bSJed Brown    test:
349c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
350c4762a1bSJed Brown 
351c4762a1bSJed Brown    test:
352c4762a1bSJed Brown       suffix: 2
353c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontr -snes_view
354c4762a1bSJed Brown       requires: !single
355c4762a1bSJed Brown 
356c4762a1bSJed Brown    test:
357c4762a1bSJed Brown       suffix: 3
358c4762a1bSJed 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
359c4762a1bSJed Brown 
36041ba4c6cSHeeho Park    test:
36141ba4c6cSHeeho Park       suffix: 4
36241ba4c6cSHeeho Park       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontrdc -snes_view
36341ba4c6cSHeeho Park       requires: !single
36441ba4c6cSHeeho Park 
365c4762a1bSJed Brown TEST*/
366