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